1454e267fSLisandro Dalcin 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3454e267fSLisandro Dalcin 4*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_1D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 5*d71ae5a4SJacob Faibussowitsch { 6454e267fSLisandro Dalcin DM_DA *da = (DM_DA *)dm->data; 7454e267fSLisandro Dalcin PetscInt i, xs, xe, Xs, Xe; 8454e267fSLisandro Dalcin PetscInt cnt = 0; 95fd66863SKarl Rupp 10454e267fSLisandro Dalcin PetscFunctionBegin; 11454e267fSLisandro Dalcin if (!da->e) { 12f951a8fcSStefano Zampini PetscInt corners[2]; 13f951a8fcSStefano Zampini 147a8be351SBarry Smith PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width"); 159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xe, NULL, NULL)); 169566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm, &Xs, NULL, NULL, &Xe, NULL, NULL)); 179371c9d4SSatish Balay xe += xs; 189371c9d4SSatish Balay Xe += Xs; 199371c9d4SSatish Balay if (xs != Xs) xs -= 1; 20454e267fSLisandro Dalcin da->ne = 1 * (xe - xs - 1); 219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + 2 * da->ne, &da->e)); 22454e267fSLisandro Dalcin for (i = xs; i < xe - 1; i++) { 23454e267fSLisandro Dalcin da->e[cnt++] = (i - Xs); 24454e267fSLisandro Dalcin da->e[cnt++] = (i - Xs + 1); 25454e267fSLisandro Dalcin } 26cfbed8a3SStefano Zampini da->nen = 2; 27cfbed8a3SStefano Zampini 28f951a8fcSStefano Zampini corners[0] = (xs - Xs); 29f951a8fcSStefano Zampini corners[1] = (xe - 1 - Xs); 309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, corners, PETSC_COPY_VALUES, &da->ecorners)); 31454e267fSLisandro Dalcin } 32454e267fSLisandro Dalcin *nel = da->ne; 33cfbed8a3SStefano Zampini *nen = da->nen; 34454e267fSLisandro Dalcin *e = da->e; 35454e267fSLisandro Dalcin PetscFunctionReturn(0); 36454e267fSLisandro Dalcin } 37454e267fSLisandro Dalcin 38*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_2D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 39*d71ae5a4SJacob Faibussowitsch { 40454e267fSLisandro Dalcin DM_DA *da = (DM_DA *)dm->data; 41454e267fSLisandro Dalcin PetscInt i, xs, xe, Xs, Xe; 42454e267fSLisandro Dalcin PetscInt j, ys, ye, Ys, Ye; 43cfbed8a3SStefano Zampini PetscInt cnt = 0, cell[4], ns = 2; 449371c9d4SSatish Balay PetscInt c, split[] = {0, 1, 3, 2, 3, 1}; 455fd66863SKarl Rupp 46454e267fSLisandro Dalcin PetscFunctionBegin; 47454e267fSLisandro Dalcin if (!da->e) { 48cfbed8a3SStefano Zampini PetscInt corners[4], nn = 0; 49f951a8fcSStefano Zampini 507a8be351SBarry Smith PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width"); 51cfbed8a3SStefano Zampini 52cfbed8a3SStefano Zampini switch (da->elementtype) { 53*d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_Q1: 54*d71ae5a4SJacob Faibussowitsch da->nen = 4; 55*d71ae5a4SJacob Faibussowitsch break; 56*d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_P1: 57*d71ae5a4SJacob Faibussowitsch da->nen = 3; 58*d71ae5a4SJacob Faibussowitsch break; 59*d71ae5a4SJacob Faibussowitsch default: 60*d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype); 61cfbed8a3SStefano Zampini } 62cfbed8a3SStefano Zampini nn = da->nen; 63cfbed8a3SStefano Zampini 64ad540459SPierre Jolivet if (da->elementtype == DMDA_ELEMENT_P1) ns = 2; 65ad540459SPierre Jolivet if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1; 669566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, &ys, NULL, &xe, &ye, NULL)); 679566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, NULL, &Xe, &Ye, NULL)); 689371c9d4SSatish Balay xe += xs; 699371c9d4SSatish Balay Xe += Xs; 709371c9d4SSatish Balay if (xs != Xs) xs -= 1; 719371c9d4SSatish Balay ye += ys; 729371c9d4SSatish Balay Ye += Ys; 739371c9d4SSatish Balay if (ys != Ys) ys -= 1; 74454e267fSLisandro Dalcin da->ne = ns * (xe - xs - 1) * (ye - ys - 1); 759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e)); 76454e267fSLisandro Dalcin for (j = ys; j < ye - 1; j++) { 77454e267fSLisandro Dalcin for (i = xs; i < xe - 1; i++) { 78454e267fSLisandro Dalcin cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs); 79454e267fSLisandro Dalcin cell[1] = (i - Xs + 1) + (j - Ys) * (Xe - Xs); 80454e267fSLisandro Dalcin cell[2] = (i - Xs + 1) + (j - Ys + 1) * (Xe - Xs); 81454e267fSLisandro Dalcin cell[3] = (i - Xs) + (j - Ys + 1) * (Xe - Xs); 82454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 838865f1eaSKarl Rupp for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]]; 84454e267fSLisandro Dalcin } 85454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 868865f1eaSKarl Rupp for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c]; 87454e267fSLisandro Dalcin } 88454e267fSLisandro Dalcin } 89454e267fSLisandro Dalcin } 90cfbed8a3SStefano Zampini 91f951a8fcSStefano Zampini corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs); 92f951a8fcSStefano Zampini corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs); 93f951a8fcSStefano Zampini corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs); 94f951a8fcSStefano Zampini corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs); 959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 4, corners, PETSC_COPY_VALUES, &da->ecorners)); 96454e267fSLisandro Dalcin } 97454e267fSLisandro Dalcin *nel = da->ne; 98cfbed8a3SStefano Zampini *nen = da->nen; 99454e267fSLisandro Dalcin *e = da->e; 100454e267fSLisandro Dalcin PetscFunctionReturn(0); 101454e267fSLisandro Dalcin } 102454e267fSLisandro Dalcin 103*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_3D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 104*d71ae5a4SJacob Faibussowitsch { 105454e267fSLisandro Dalcin DM_DA *da = (DM_DA *)dm->data; 106454e267fSLisandro Dalcin PetscInt i, xs, xe, Xs, Xe; 107454e267fSLisandro Dalcin PetscInt j, ys, ye, Ys, Ye; 108454e267fSLisandro Dalcin PetscInt k, zs, ze, Zs, Ze; 109cfbed8a3SStefano Zampini PetscInt cnt = 0, cell[8], ns = 6; 1109371c9d4SSatish Balay PetscInt c, split[] = {0, 1, 3, 7, 0, 1, 7, 4, 1, 2, 3, 7, 1, 2, 7, 6, 1, 4, 5, 7, 1, 5, 6, 7}; 1115fd66863SKarl Rupp 112454e267fSLisandro Dalcin PetscFunctionBegin; 113454e267fSLisandro Dalcin if (!da->e) { 114cfbed8a3SStefano Zampini PetscInt corners[8], nn = 0; 115f951a8fcSStefano Zampini 1167a8be351SBarry Smith PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width"); 117cfbed8a3SStefano Zampini 118cfbed8a3SStefano Zampini switch (da->elementtype) { 119*d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_Q1: 120*d71ae5a4SJacob Faibussowitsch da->nen = 8; 121*d71ae5a4SJacob Faibussowitsch break; 122*d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_P1: 123*d71ae5a4SJacob Faibussowitsch da->nen = 4; 124*d71ae5a4SJacob Faibussowitsch break; 125*d71ae5a4SJacob Faibussowitsch default: 126*d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype); 127cfbed8a3SStefano Zampini } 128cfbed8a3SStefano Zampini nn = da->nen; 129cfbed8a3SStefano Zampini 130ad540459SPierre Jolivet if (da->elementtype == DMDA_ELEMENT_P1) ns = 6; 131ad540459SPierre Jolivet if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1; 1329566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, &ys, &zs, &xe, &ye, &ze)); 1339566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze)); 1349371c9d4SSatish Balay xe += xs; 1359371c9d4SSatish Balay Xe += Xs; 1369371c9d4SSatish Balay if (xs != Xs) xs -= 1; 1379371c9d4SSatish Balay ye += ys; 1389371c9d4SSatish Balay Ye += Ys; 1399371c9d4SSatish Balay if (ys != Ys) ys -= 1; 1409371c9d4SSatish Balay ze += zs; 1419371c9d4SSatish Balay Ze += Zs; 1429371c9d4SSatish Balay if (zs != Zs) zs -= 1; 143454e267fSLisandro Dalcin da->ne = ns * (xe - xs - 1) * (ye - ys - 1) * (ze - zs - 1); 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e)); 145454e267fSLisandro Dalcin for (k = zs; k < ze - 1; k++) { 146454e267fSLisandro Dalcin for (j = ys; j < ye - 1; j++) { 147454e267fSLisandro Dalcin for (i = xs; i < xe - 1; i++) { 148454e267fSLisandro Dalcin cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys); 1492da392ccSBarry Smith cell[1] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys); 1502da392ccSBarry Smith cell[2] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys); 1512da392ccSBarry Smith cell[3] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys); 1522da392ccSBarry Smith cell[4] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys); 1532da392ccSBarry Smith cell[5] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys); 1542da392ccSBarry Smith cell[6] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys); 1552da392ccSBarry Smith cell[7] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys); 156454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 1578865f1eaSKarl Rupp for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]]; 158454e267fSLisandro Dalcin } 159454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 1608865f1eaSKarl Rupp for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c]; 161454e267fSLisandro Dalcin } 162454e267fSLisandro Dalcin } 163454e267fSLisandro Dalcin } 164454e267fSLisandro Dalcin } 165cfbed8a3SStefano Zampini 166f951a8fcSStefano Zampini corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 167f951a8fcSStefano Zampini corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 168f951a8fcSStefano Zampini corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 169f951a8fcSStefano Zampini corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 170f951a8fcSStefano Zampini corners[4] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys); 171f951a8fcSStefano Zampini corners[5] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys); 172f951a8fcSStefano Zampini corners[6] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys); 173f951a8fcSStefano Zampini corners[7] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys); 1749566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 8, corners, PETSC_COPY_VALUES, &da->ecorners)); 175454e267fSLisandro Dalcin } 176454e267fSLisandro Dalcin *nel = da->ne; 177cfbed8a3SStefano Zampini *nen = da->nen; 178454e267fSLisandro Dalcin *e = da->e; 179454e267fSLisandro Dalcin PetscFunctionReturn(0); 180454e267fSLisandro Dalcin } 181454e267fSLisandro Dalcin 182f5f57ec0SBarry Smith /*@ 183d4a6ed37SStefano Zampini DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left 184d4a6ed37SStefano Zampini corner of the non-overlapping decomposition identified by DMDAGetElements() 185d4a6ed37SStefano Zampini 186d4a6ed37SStefano Zampini Not Collective 187d4a6ed37SStefano Zampini 188d4a6ed37SStefano Zampini Input Parameter: 189d4a6ed37SStefano Zampini . da - the DM object 190d4a6ed37SStefano Zampini 191d4a6ed37SStefano Zampini Output Parameters: 192d4a6ed37SStefano Zampini + gx - the x index 193d4a6ed37SStefano Zampini . gy - the y index 194d4a6ed37SStefano Zampini - gz - the z index 195d4a6ed37SStefano Zampini 196d4a6ed37SStefano Zampini Level: intermediate 197d4a6ed37SStefano Zampini 198d4a6ed37SStefano Zampini Notes: 199d4a6ed37SStefano Zampini 200db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()` 201d4a6ed37SStefano Zampini @*/ 202*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz) 203*d71ae5a4SJacob Faibussowitsch { 204d4a6ed37SStefano Zampini PetscInt xs, Xs; 205d4a6ed37SStefano Zampini PetscInt ys, Ys; 206d4a6ed37SStefano Zampini PetscInt zs, Zs; 207d4a6ed37SStefano Zampini PetscBool isda; 208d4a6ed37SStefano Zampini 209d4a6ed37SStefano Zampini PetscFunctionBegin; 210a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 211d4a6ed37SStefano Zampini if (gx) PetscValidIntPointer(gx, 2); 212d4a6ed37SStefano Zampini if (gy) PetscValidIntPointer(gy, 3); 213d4a6ed37SStefano Zampini if (gz) PetscValidIntPointer(gz, 4); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 2157a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 2169566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, NULL, NULL, NULL)); 2179566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL)); 218d4a6ed37SStefano Zampini if (xs != Xs) xs -= 1; 219d4a6ed37SStefano Zampini if (ys != Ys) ys -= 1; 220d4a6ed37SStefano Zampini if (zs != Zs) zs -= 1; 221d4a6ed37SStefano Zampini if (gx) *gx = xs; 222d4a6ed37SStefano Zampini if (gy) *gy = ys; 223d4a6ed37SStefano Zampini if (gz) *gz = zs; 224d4a6ed37SStefano Zampini PetscFunctionReturn(0); 225d4a6ed37SStefano Zampini } 226d4a6ed37SStefano Zampini 227d4a6ed37SStefano Zampini /*@ 228d4a6ed37SStefano Zampini DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by DMDAGetElements() 229d4a6ed37SStefano Zampini 230d4a6ed37SStefano Zampini Not Collective 231d4a6ed37SStefano Zampini 232d4a6ed37SStefano Zampini Input Parameter: 233d4a6ed37SStefano Zampini . da - the DM object 234d4a6ed37SStefano Zampini 235d4a6ed37SStefano Zampini Output Parameters: 236d4a6ed37SStefano Zampini + mx - number of local elements in x-direction 237d4a6ed37SStefano Zampini . my - number of local elements in y-direction 238d4a6ed37SStefano Zampini - mz - number of local elements in z-direction 239d4a6ed37SStefano Zampini 240d4a6ed37SStefano Zampini Level: intermediate 241d4a6ed37SStefano Zampini 24295452b02SPatrick Sanan Notes: 24395452b02SPatrick Sanan It returns the same number of elements, irrespective of the DMDAElementType 244d4a6ed37SStefano Zampini 245db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements` 246d4a6ed37SStefano Zampini @*/ 247*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz) 248*d71ae5a4SJacob Faibussowitsch { 249d4a6ed37SStefano Zampini PetscInt xs, xe, Xs; 250d4a6ed37SStefano Zampini PetscInt ys, ye, Ys; 251d4a6ed37SStefano Zampini PetscInt zs, ze, Zs; 252d4a6ed37SStefano Zampini PetscInt dim; 253d4a6ed37SStefano Zampini PetscBool isda; 254d4a6ed37SStefano Zampini 255d4a6ed37SStefano Zampini PetscFunctionBegin; 256a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 257d4a6ed37SStefano Zampini if (mx) PetscValidIntPointer(mx, 2); 258d4a6ed37SStefano Zampini if (my) PetscValidIntPointer(my, 3); 259d4a6ed37SStefano Zampini if (mz) PetscValidIntPointer(mz, 4); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 2617a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 2629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xe, &ye, &ze)); 2639566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL)); 2649371c9d4SSatish Balay xe += xs; 2659371c9d4SSatish Balay if (xs != Xs) xs -= 1; 2669371c9d4SSatish Balay ye += ys; 2679371c9d4SSatish Balay if (ys != Ys) ys -= 1; 2689371c9d4SSatish Balay ze += zs; 2699371c9d4SSatish Balay if (zs != Zs) zs -= 1; 270d4a6ed37SStefano Zampini if (mx) *mx = 0; 271d4a6ed37SStefano Zampini if (my) *my = 0; 272d4a6ed37SStefano Zampini if (mz) *mz = 0; 2739566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 274d4a6ed37SStefano Zampini switch (dim) { 275d4a6ed37SStefano Zampini case 3: 276f4d061e9SPierre Jolivet if (mz) *mz = ze - zs - 1; /* fall through */ 277d4a6ed37SStefano Zampini case 2: 278f4d061e9SPierre Jolivet if (my) *my = ye - ys - 1; /* fall through */ 279d4a6ed37SStefano Zampini case 1: 280d4a6ed37SStefano Zampini if (mx) *mx = xe - xs - 1; 281d4a6ed37SStefano Zampini break; 282d4a6ed37SStefano Zampini } 283d4a6ed37SStefano Zampini PetscFunctionReturn(0); 284d4a6ed37SStefano Zampini } 285d4a6ed37SStefano Zampini 286d4a6ed37SStefano Zampini /*@ 2872dde6fd4SLisandro Dalcin DMDASetElementType - Sets the element type to be returned by DMDAGetElements() 2882dde6fd4SLisandro Dalcin 2892dde6fd4SLisandro Dalcin Not Collective 2902dde6fd4SLisandro Dalcin 2912dde6fd4SLisandro Dalcin Input Parameter: 2922dde6fd4SLisandro Dalcin . da - the DMDA object 2932dde6fd4SLisandro Dalcin 2942dde6fd4SLisandro Dalcin Output Parameters: 2956420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 2962dde6fd4SLisandro Dalcin 2972dde6fd4SLisandro Dalcin Level: intermediate 2982dde6fd4SLisandro Dalcin 299db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDAGetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()` 3002dde6fd4SLisandro Dalcin @*/ 301*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype) 302*d71ae5a4SJacob Faibussowitsch { 3032dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 304f951a8fcSStefano Zampini PetscBool isda; 3052dde6fd4SLisandro Dalcin 3062dde6fd4SLisandro Dalcin PetscFunctionBegin; 307a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 3082dde6fd4SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, etype, 2); 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 310f951a8fcSStefano Zampini if (!isda) PetscFunctionReturn(0); 3112dde6fd4SLisandro Dalcin if (dd->elementtype != etype) { 3129566063dSJacob Faibussowitsch PetscCall(PetscFree(dd->e)); 3139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dd->ecorners)); 3148865f1eaSKarl Rupp 3152dde6fd4SLisandro Dalcin dd->elementtype = etype; 3162dde6fd4SLisandro Dalcin dd->ne = 0; 317cfbed8a3SStefano Zampini dd->nen = 0; 3180298fd71SBarry Smith dd->e = NULL; 3192dde6fd4SLisandro Dalcin } 3202dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3212dde6fd4SLisandro Dalcin } 3222dde6fd4SLisandro Dalcin 323f5f57ec0SBarry Smith /*@ 3242dde6fd4SLisandro Dalcin DMDAGetElementType - Gets the element type to be returned by DMDAGetElements() 3252dde6fd4SLisandro Dalcin 3262dde6fd4SLisandro Dalcin Not Collective 3272dde6fd4SLisandro Dalcin 3282dde6fd4SLisandro Dalcin Input Parameter: 3292dde6fd4SLisandro Dalcin . da - the DMDA object 3302dde6fd4SLisandro Dalcin 3312dde6fd4SLisandro Dalcin Output Parameters: 3326420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 3332dde6fd4SLisandro Dalcin 3342dde6fd4SLisandro Dalcin Level: intermediate 3352dde6fd4SLisandro Dalcin 336db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()` 3372dde6fd4SLisandro Dalcin @*/ 338*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype) 339*d71ae5a4SJacob Faibussowitsch { 3402dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 341f951a8fcSStefano Zampini PetscBool isda; 3425fd66863SKarl Rupp 3432dde6fd4SLisandro Dalcin PetscFunctionBegin; 344a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 3452dde6fd4SLisandro Dalcin PetscValidPointer(etype, 2); 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 3477a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 3482dde6fd4SLisandro Dalcin *etype = dd->elementtype; 3492dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3502dde6fd4SLisandro Dalcin } 3512dde6fd4SLisandro Dalcin 3522dde6fd4SLisandro Dalcin /*@C 3532dde6fd4SLisandro Dalcin DMDAGetElements - Gets an array containing the indices (in local coordinates) 3542dde6fd4SLisandro Dalcin of all the local elements 3552dde6fd4SLisandro Dalcin 3562dde6fd4SLisandro Dalcin Not Collective 3572dde6fd4SLisandro Dalcin 3582dde6fd4SLisandro Dalcin Input Parameter: 3592dde6fd4SLisandro Dalcin . dm - the DM object 3602dde6fd4SLisandro Dalcin 3612dde6fd4SLisandro Dalcin Output Parameters: 3622dde6fd4SLisandro Dalcin + nel - number of local elements 3632dde6fd4SLisandro Dalcin . nen - number of element nodes 364c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 3652dde6fd4SLisandro Dalcin 3662dde6fd4SLisandro Dalcin Level: intermediate 3672dde6fd4SLisandro Dalcin 36893386343SJed Brown Notes: 36993386343SJed Brown Call DMDARestoreElements() once you have finished accessing the elements. 37093386343SJed Brown 37145950774SSatish Balay Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes. 372009f0576SBarry Smith 373009f0576SBarry Smith If on each process you integrate over its owned elements and use ADD_VALUES in Vec/MatSetValuesLocal() then you'll obtain the correct result. 374009f0576SBarry Smith 375f5f57ec0SBarry Smith Not supported in Fortran 376f5f57ec0SBarry Smith 377db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`, `DMGlobalToLocalBegin()`, `DMLocalToGlobalBegin()` 3782dde6fd4SLisandro Dalcin @*/ 379*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 380*d71ae5a4SJacob Faibussowitsch { 381c73cfb54SMatthew G. Knepley PetscInt dim; 38299c57625SBarry Smith DM_DA *dd = (DM_DA *)dm->data; 383f951a8fcSStefano Zampini PetscBool isda; 3845fd66863SKarl Rupp 385454e267fSLisandro Dalcin PetscFunctionBegin; 386a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 387f951a8fcSStefano Zampini PetscValidIntPointer(nel, 2); 388f951a8fcSStefano Zampini PetscValidIntPointer(nen, 3); 389f951a8fcSStefano Zampini PetscValidPointer(e, 4); 3909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 3917a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name); 39208401ef6SPierre Jolivet PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMDAGetElements() requires you use a stencil type of DMDA_STENCIL_BOX"); 3939566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 394cfbed8a3SStefano Zampini if (dd->e) { 395cfbed8a3SStefano Zampini *nel = dd->ne; 396cfbed8a3SStefano Zampini *nen = dd->nen; 397cfbed8a3SStefano Zampini *e = dd->e; 398cfbed8a3SStefano Zampini PetscFunctionReturn(0); 399cfbed8a3SStefano Zampini } 400c73cfb54SMatthew G. Knepley if (dim == -1) { 4019371c9d4SSatish Balay *nel = 0; 4029371c9d4SSatish Balay *nen = 0; 4039371c9d4SSatish Balay *e = NULL; 404c73cfb54SMatthew G. Knepley } else if (dim == 1) { 4059566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_1D(dm, nel, nen, e)); 406c73cfb54SMatthew G. Knepley } else if (dim == 2) { 4079566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_2D(dm, nel, nen, e)); 408c73cfb54SMatthew G. Knepley } else if (dim == 3) { 4099566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_3D(dm, nel, nen, e)); 41063a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 411454e267fSLisandro Dalcin PetscFunctionReturn(0); 412454e267fSLisandro Dalcin } 4132dde6fd4SLisandro Dalcin 414f951a8fcSStefano Zampini /*@ 415d4a6ed37SStefano Zampini DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates) 416f951a8fcSStefano Zampini of the non-overlapping decomposition identified by DMDAGetElements 417f951a8fcSStefano Zampini 418f951a8fcSStefano Zampini Not Collective 419f951a8fcSStefano Zampini 420f951a8fcSStefano Zampini Input Parameter: 421f951a8fcSStefano Zampini . dm - the DM object 422f951a8fcSStefano Zampini 423f951a8fcSStefano Zampini Output Parameters: 424f951a8fcSStefano Zampini . is - the index set 425f951a8fcSStefano Zampini 426f951a8fcSStefano Zampini Level: intermediate 427f951a8fcSStefano Zampini 42895452b02SPatrick Sanan Notes: 42995452b02SPatrick Sanan Call DMDARestoreSubdomainCornersIS() once you have finished accessing the index set. 430f951a8fcSStefano Zampini 431db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElementsCornersIS()` 432f951a8fcSStefano Zampini @*/ 433*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetSubdomainCornersIS(DM dm, IS *is) 434*d71ae5a4SJacob Faibussowitsch { 435f951a8fcSStefano Zampini DM_DA *dd = (DM_DA *)dm->data; 436f951a8fcSStefano Zampini PetscBool isda; 437f951a8fcSStefano Zampini 438f951a8fcSStefano Zampini PetscFunctionBegin; 439a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 440f951a8fcSStefano Zampini PetscValidPointer(is, 2); 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 4427a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name); 44308401ef6SPierre Jolivet PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMDAGetElement() requires you use a stencil type of DMDA_STENCIL_BOX"); 444f951a8fcSStefano Zampini if (!dd->ecorners) { /* compute elements if not yet done */ 445f951a8fcSStefano Zampini const PetscInt *e; 446f951a8fcSStefano Zampini PetscInt nel, nen; 447f951a8fcSStefano Zampini 4489566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dm, &nel, &nen, &e)); 4499566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dm, &nel, &nen, &e)); 450f951a8fcSStefano Zampini } 451f951a8fcSStefano Zampini *is = dd->ecorners; 452f951a8fcSStefano Zampini PetscFunctionReturn(0); 453f951a8fcSStefano Zampini } 454f951a8fcSStefano Zampini 4552dde6fd4SLisandro Dalcin /*@C 456009f0576SBarry Smith DMDARestoreElements - Restores the array obtained with DMDAGetElements() 4572dde6fd4SLisandro Dalcin 4582dde6fd4SLisandro Dalcin Not Collective 4592dde6fd4SLisandro Dalcin 460d8d19677SJose E. Roman Input Parameters: 4612dde6fd4SLisandro Dalcin + dm - the DM object 4622dde6fd4SLisandro Dalcin . nel - number of local elements 4632dde6fd4SLisandro Dalcin . nen - number of element nodes 464c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 4652dde6fd4SLisandro Dalcin 4662dde6fd4SLisandro Dalcin Level: intermediate 4672dde6fd4SLisandro Dalcin 468009f0576SBarry Smith Note: You should not access these values after you have called this routine. 469009f0576SBarry Smith 470c67aec63SBarry Smith This restore signals the DMDA object that you no longer need access to the array information. 471009f0576SBarry Smith 472f5f57ec0SBarry Smith Not supported in Fortran 473f5f57ec0SBarry Smith 474db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()` 4752dde6fd4SLisandro Dalcin @*/ 476*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 477*d71ae5a4SJacob Faibussowitsch { 4782dde6fd4SLisandro Dalcin PetscFunctionBegin; 479a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 4802dde6fd4SLisandro Dalcin PetscValidIntPointer(nel, 2); 4812dde6fd4SLisandro Dalcin PetscValidIntPointer(nen, 3); 4822dde6fd4SLisandro Dalcin PetscValidPointer(e, 4); 4839cc8bc54SJed Brown *nel = 0; 4849cc8bc54SJed Brown *nen = -1; 4859cc8bc54SJed Brown *e = NULL; 4862dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 4872dde6fd4SLisandro Dalcin } 488f951a8fcSStefano Zampini 489f951a8fcSStefano Zampini /*@ 490d4a6ed37SStefano Zampini DMDARestoreSubdomainCornersIS - Restores the IS obtained with DMDAGetSubdomainCornersIS() 491f951a8fcSStefano Zampini 492f951a8fcSStefano Zampini Not Collective 493f951a8fcSStefano Zampini 494d8d19677SJose E. Roman Input Parameters: 495f951a8fcSStefano Zampini + dm - the DM object 496f951a8fcSStefano Zampini - is - the index set 497f951a8fcSStefano Zampini 498f951a8fcSStefano Zampini Level: intermediate 499f951a8fcSStefano Zampini 500f951a8fcSStefano Zampini Note: 501f951a8fcSStefano Zampini 502db781477SPatrick Sanan .seealso: `DMDAElementType`, `DMDASetElementType()`, `DMDAGetSubdomainCornersIS()` 503f951a8fcSStefano Zampini @*/ 504*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm, IS *is) 505*d71ae5a4SJacob Faibussowitsch { 506f951a8fcSStefano Zampini PetscFunctionBegin; 507a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 508f951a8fcSStefano Zampini PetscValidHeaderSpecific(*is, IS_CLASSID, 2); 509f951a8fcSStefano Zampini *is = NULL; 510f951a8fcSStefano Zampini PetscFunctionReturn(0); 511f951a8fcSStefano Zampini } 512