1454e267fSLisandro Dalcin 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3454e267fSLisandro Dalcin 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_1D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 5d71ae5a4SJacob 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 38d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_2D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 39d71ae5a4SJacob 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) { 53d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_Q1: 54d71ae5a4SJacob Faibussowitsch da->nen = 4; 55d71ae5a4SJacob Faibussowitsch break; 56d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_P1: 57d71ae5a4SJacob Faibussowitsch da->nen = 3; 58d71ae5a4SJacob Faibussowitsch break; 59d71ae5a4SJacob Faibussowitsch default: 60d71ae5a4SJacob 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 103d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_3D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 104d71ae5a4SJacob 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) { 119d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_Q1: 120d71ae5a4SJacob Faibussowitsch da->nen = 8; 121d71ae5a4SJacob Faibussowitsch break; 122d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_P1: 123d71ae5a4SJacob Faibussowitsch da->nen = 4; 124d71ae5a4SJacob Faibussowitsch break; 125d71ae5a4SJacob Faibussowitsch default: 126d71ae5a4SJacob 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 184*dce8aebaSBarry Smith corner of the non-overlapping decomposition identified by `DMDAGetElements()` 185d4a6ed37SStefano Zampini 186d4a6ed37SStefano Zampini Not Collective 187d4a6ed37SStefano Zampini 188d4a6ed37SStefano Zampini Input Parameter: 189*dce8aebaSBarry Smith . da - the `DMDA` 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 198*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()` 199d4a6ed37SStefano Zampini @*/ 200d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz) 201d71ae5a4SJacob Faibussowitsch { 202d4a6ed37SStefano Zampini PetscInt xs, Xs; 203d4a6ed37SStefano Zampini PetscInt ys, Ys; 204d4a6ed37SStefano Zampini PetscInt zs, Zs; 205d4a6ed37SStefano Zampini PetscBool isda; 206d4a6ed37SStefano Zampini 207d4a6ed37SStefano Zampini PetscFunctionBegin; 208a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 209d4a6ed37SStefano Zampini if (gx) PetscValidIntPointer(gx, 2); 210d4a6ed37SStefano Zampini if (gy) PetscValidIntPointer(gy, 3); 211d4a6ed37SStefano Zampini if (gz) PetscValidIntPointer(gz, 4); 2129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 2137a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 2149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, NULL, NULL, NULL)); 2159566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL)); 216d4a6ed37SStefano Zampini if (xs != Xs) xs -= 1; 217d4a6ed37SStefano Zampini if (ys != Ys) ys -= 1; 218d4a6ed37SStefano Zampini if (zs != Zs) zs -= 1; 219d4a6ed37SStefano Zampini if (gx) *gx = xs; 220d4a6ed37SStefano Zampini if (gy) *gy = ys; 221d4a6ed37SStefano Zampini if (gz) *gz = zs; 222d4a6ed37SStefano Zampini PetscFunctionReturn(0); 223d4a6ed37SStefano Zampini } 224d4a6ed37SStefano Zampini 225d4a6ed37SStefano Zampini /*@ 226*dce8aebaSBarry Smith DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by `DMDAGetElements()` 227d4a6ed37SStefano Zampini 228d4a6ed37SStefano Zampini Not Collective 229d4a6ed37SStefano Zampini 230d4a6ed37SStefano Zampini Input Parameter: 231*dce8aebaSBarry Smith . da - the `DMDA` object 232d4a6ed37SStefano Zampini 233d4a6ed37SStefano Zampini Output Parameters: 234d4a6ed37SStefano Zampini + mx - number of local elements in x-direction 235d4a6ed37SStefano Zampini . my - number of local elements in y-direction 236d4a6ed37SStefano Zampini - mz - number of local elements in z-direction 237d4a6ed37SStefano Zampini 238d4a6ed37SStefano Zampini Level: intermediate 239d4a6ed37SStefano Zampini 240*dce8aebaSBarry Smith Note: 241*dce8aebaSBarry Smith It returns the same number of elements, irrespective of the `DMDAElementType` 242d4a6ed37SStefano Zampini 243*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements` 244d4a6ed37SStefano Zampini @*/ 245d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz) 246d71ae5a4SJacob Faibussowitsch { 247d4a6ed37SStefano Zampini PetscInt xs, xe, Xs; 248d4a6ed37SStefano Zampini PetscInt ys, ye, Ys; 249d4a6ed37SStefano Zampini PetscInt zs, ze, Zs; 250d4a6ed37SStefano Zampini PetscInt dim; 251d4a6ed37SStefano Zampini PetscBool isda; 252d4a6ed37SStefano Zampini 253d4a6ed37SStefano Zampini PetscFunctionBegin; 254a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 255d4a6ed37SStefano Zampini if (mx) PetscValidIntPointer(mx, 2); 256d4a6ed37SStefano Zampini if (my) PetscValidIntPointer(my, 3); 257d4a6ed37SStefano Zampini if (mz) PetscValidIntPointer(mz, 4); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 2597a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 2609566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xe, &ye, &ze)); 2619566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, NULL, NULL, NULL)); 2629371c9d4SSatish Balay xe += xs; 2639371c9d4SSatish Balay if (xs != Xs) xs -= 1; 2649371c9d4SSatish Balay ye += ys; 2659371c9d4SSatish Balay if (ys != Ys) ys -= 1; 2669371c9d4SSatish Balay ze += zs; 2679371c9d4SSatish Balay if (zs != Zs) zs -= 1; 268d4a6ed37SStefano Zampini if (mx) *mx = 0; 269d4a6ed37SStefano Zampini if (my) *my = 0; 270d4a6ed37SStefano Zampini if (mz) *mz = 0; 2719566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 272d4a6ed37SStefano Zampini switch (dim) { 273d4a6ed37SStefano Zampini case 3: 274f4d061e9SPierre Jolivet if (mz) *mz = ze - zs - 1; /* fall through */ 275d4a6ed37SStefano Zampini case 2: 276f4d061e9SPierre Jolivet if (my) *my = ye - ys - 1; /* fall through */ 277d4a6ed37SStefano Zampini case 1: 278d4a6ed37SStefano Zampini if (mx) *mx = xe - xs - 1; 279d4a6ed37SStefano Zampini break; 280d4a6ed37SStefano Zampini } 281d4a6ed37SStefano Zampini PetscFunctionReturn(0); 282d4a6ed37SStefano Zampini } 283d4a6ed37SStefano Zampini 284d4a6ed37SStefano Zampini /*@ 285*dce8aebaSBarry Smith DMDASetElementType - Sets the element type to be returned by `DMDAGetElements()` 2862dde6fd4SLisandro Dalcin 2872dde6fd4SLisandro Dalcin Not Collective 2882dde6fd4SLisandro Dalcin 2892dde6fd4SLisandro Dalcin Input Parameter: 290*dce8aebaSBarry Smith . da - the `DMDA` object 2912dde6fd4SLisandro Dalcin 2922dde6fd4SLisandro Dalcin Output Parameters: 293*dce8aebaSBarry Smith . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1` 2942dde6fd4SLisandro Dalcin 2952dde6fd4SLisandro Dalcin Level: intermediate 2962dde6fd4SLisandro Dalcin 297*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDAGetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()` 2982dde6fd4SLisandro Dalcin @*/ 299d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype) 300d71ae5a4SJacob Faibussowitsch { 3012dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 302f951a8fcSStefano Zampini PetscBool isda; 3032dde6fd4SLisandro Dalcin 3042dde6fd4SLisandro Dalcin PetscFunctionBegin; 305a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 3062dde6fd4SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, etype, 2); 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 308f951a8fcSStefano Zampini if (!isda) PetscFunctionReturn(0); 3092dde6fd4SLisandro Dalcin if (dd->elementtype != etype) { 3109566063dSJacob Faibussowitsch PetscCall(PetscFree(dd->e)); 3119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dd->ecorners)); 3128865f1eaSKarl Rupp 3132dde6fd4SLisandro Dalcin dd->elementtype = etype; 3142dde6fd4SLisandro Dalcin dd->ne = 0; 315cfbed8a3SStefano Zampini dd->nen = 0; 3160298fd71SBarry Smith dd->e = NULL; 3172dde6fd4SLisandro Dalcin } 3182dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3192dde6fd4SLisandro Dalcin } 3202dde6fd4SLisandro Dalcin 321f5f57ec0SBarry Smith /*@ 322*dce8aebaSBarry Smith DMDAGetElementType - Gets the element type to be returned by `DMDAGetElements()` 3232dde6fd4SLisandro Dalcin 3242dde6fd4SLisandro Dalcin Not Collective 3252dde6fd4SLisandro Dalcin 3262dde6fd4SLisandro Dalcin Input Parameter: 327*dce8aebaSBarry Smith . da - the `DMDA` object 3282dde6fd4SLisandro Dalcin 3292dde6fd4SLisandro Dalcin Output Parameters: 330*dce8aebaSBarry Smith . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1` 3312dde6fd4SLisandro Dalcin 3322dde6fd4SLisandro Dalcin Level: intermediate 3332dde6fd4SLisandro Dalcin 334*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()` 3352dde6fd4SLisandro Dalcin @*/ 336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype) 337d71ae5a4SJacob Faibussowitsch { 3382dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 339f951a8fcSStefano Zampini PetscBool isda; 3405fd66863SKarl Rupp 3412dde6fd4SLisandro Dalcin PetscFunctionBegin; 342a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 3432dde6fd4SLisandro Dalcin PetscValidPointer(etype, 2); 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda)); 3457a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)da)->type_name); 3462dde6fd4SLisandro Dalcin *etype = dd->elementtype; 3472dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3482dde6fd4SLisandro Dalcin } 3492dde6fd4SLisandro Dalcin 3502dde6fd4SLisandro Dalcin /*@C 3512dde6fd4SLisandro Dalcin DMDAGetElements - Gets an array containing the indices (in local coordinates) 3522dde6fd4SLisandro Dalcin of all the local elements 3532dde6fd4SLisandro Dalcin 3542dde6fd4SLisandro Dalcin Not Collective 3552dde6fd4SLisandro Dalcin 3562dde6fd4SLisandro Dalcin Input Parameter: 357*dce8aebaSBarry Smith . dm - the `DMDA` object 3582dde6fd4SLisandro Dalcin 3592dde6fd4SLisandro Dalcin Output Parameters: 3602dde6fd4SLisandro Dalcin + nel - number of local elements 3612dde6fd4SLisandro Dalcin . nen - number of element nodes 362c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 3632dde6fd4SLisandro Dalcin 3642dde6fd4SLisandro Dalcin Level: intermediate 3652dde6fd4SLisandro Dalcin 36693386343SJed Brown Notes: 367*dce8aebaSBarry Smith Call `DMDARestoreElements()` once you have finished accessing the elements. 36893386343SJed Brown 36945950774SSatish Balay Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes. 370009f0576SBarry Smith 371*dce8aebaSBarry 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. 372009f0576SBarry Smith 373*dce8aebaSBarry Smith Fortran Note: 374f5f57ec0SBarry Smith Not supported in Fortran 375f5f57ec0SBarry Smith 376*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`, `DMGlobalToLocalBegin()`, `DMLocalToGlobalBegin()` 3772dde6fd4SLisandro Dalcin @*/ 378d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 379d71ae5a4SJacob Faibussowitsch { 380c73cfb54SMatthew G. Knepley PetscInt dim; 38199c57625SBarry Smith DM_DA *dd = (DM_DA *)dm->data; 382f951a8fcSStefano Zampini PetscBool isda; 3835fd66863SKarl Rupp 384454e267fSLisandro Dalcin PetscFunctionBegin; 385a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 386f951a8fcSStefano Zampini PetscValidIntPointer(nel, 2); 387f951a8fcSStefano Zampini PetscValidIntPointer(nen, 3); 388f951a8fcSStefano Zampini PetscValidPointer(e, 4); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 3907a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name); 39108401ef6SPierre 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"); 3929566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 393cfbed8a3SStefano Zampini if (dd->e) { 394cfbed8a3SStefano Zampini *nel = dd->ne; 395cfbed8a3SStefano Zampini *nen = dd->nen; 396cfbed8a3SStefano Zampini *e = dd->e; 397cfbed8a3SStefano Zampini PetscFunctionReturn(0); 398cfbed8a3SStefano Zampini } 399c73cfb54SMatthew G. Knepley if (dim == -1) { 4009371c9d4SSatish Balay *nel = 0; 4019371c9d4SSatish Balay *nen = 0; 4029371c9d4SSatish Balay *e = NULL; 403c73cfb54SMatthew G. Knepley } else if (dim == 1) { 4049566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_1D(dm, nel, nen, e)); 405c73cfb54SMatthew G. Knepley } else if (dim == 2) { 4069566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_2D(dm, nel, nen, e)); 407c73cfb54SMatthew G. Knepley } else if (dim == 3) { 4089566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_3D(dm, nel, nen, e)); 40963a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 410454e267fSLisandro Dalcin PetscFunctionReturn(0); 411454e267fSLisandro Dalcin } 4122dde6fd4SLisandro Dalcin 413f951a8fcSStefano Zampini /*@ 414d4a6ed37SStefano Zampini DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates) 415*dce8aebaSBarry Smith of the non-overlapping decomposition identified by `DMDAGetElements()` 416f951a8fcSStefano Zampini 417f951a8fcSStefano Zampini Not Collective 418f951a8fcSStefano Zampini 419f951a8fcSStefano Zampini Input Parameter: 420*dce8aebaSBarry Smith . dm - the `DMDA` object 421f951a8fcSStefano Zampini 422f951a8fcSStefano Zampini Output Parameters: 423f951a8fcSStefano Zampini . is - the index set 424f951a8fcSStefano Zampini 425f951a8fcSStefano Zampini Level: intermediate 426f951a8fcSStefano Zampini 427*dce8aebaSBarry Smith Note: 428*dce8aebaSBarry Smith Call `DMDARestoreSubdomainCornersIS()` once you have finished accessing the index set. 429f951a8fcSStefano Zampini 430*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElementsCornersIS()` 431f951a8fcSStefano Zampini @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetSubdomainCornersIS(DM dm, IS *is) 433d71ae5a4SJacob Faibussowitsch { 434f951a8fcSStefano Zampini DM_DA *dd = (DM_DA *)dm->data; 435f951a8fcSStefano Zampini PetscBool isda; 436f951a8fcSStefano Zampini 437f951a8fcSStefano Zampini PetscFunctionBegin; 438a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 439f951a8fcSStefano Zampini PetscValidPointer(is, 2); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 4417a8be351SBarry Smith PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name); 44208401ef6SPierre 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"); 443f951a8fcSStefano Zampini if (!dd->ecorners) { /* compute elements if not yet done */ 444f951a8fcSStefano Zampini const PetscInt *e; 445f951a8fcSStefano Zampini PetscInt nel, nen; 446f951a8fcSStefano Zampini 4479566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dm, &nel, &nen, &e)); 4489566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dm, &nel, &nen, &e)); 449f951a8fcSStefano Zampini } 450f951a8fcSStefano Zampini *is = dd->ecorners; 451f951a8fcSStefano Zampini PetscFunctionReturn(0); 452f951a8fcSStefano Zampini } 453f951a8fcSStefano Zampini 4542dde6fd4SLisandro Dalcin /*@C 455*dce8aebaSBarry Smith DMDARestoreElements - Restores the array obtained with `DMDAGetElements()` 4562dde6fd4SLisandro Dalcin 4572dde6fd4SLisandro Dalcin Not Collective 4582dde6fd4SLisandro Dalcin 459d8d19677SJose E. Roman Input Parameters: 4602dde6fd4SLisandro Dalcin + dm - the DM object 4612dde6fd4SLisandro Dalcin . nel - number of local elements 4622dde6fd4SLisandro Dalcin . nen - number of element nodes 463c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 4642dde6fd4SLisandro Dalcin 4652dde6fd4SLisandro Dalcin Level: intermediate 4662dde6fd4SLisandro Dalcin 467*dce8aebaSBarry Smith Note: 468*dce8aebaSBarry Smith This restore signals the `DMDA` object that you no longer need access to the array information. 469009f0576SBarry Smith 470*dce8aebaSBarry Smith Fortran Note: 471f5f57ec0SBarry Smith Not supported in Fortran 472f5f57ec0SBarry Smith 473*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()` 4742dde6fd4SLisandro Dalcin @*/ 475d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[]) 476d71ae5a4SJacob Faibussowitsch { 4772dde6fd4SLisandro Dalcin PetscFunctionBegin; 478a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 4792dde6fd4SLisandro Dalcin PetscValidIntPointer(nel, 2); 4802dde6fd4SLisandro Dalcin PetscValidIntPointer(nen, 3); 4812dde6fd4SLisandro Dalcin PetscValidPointer(e, 4); 4829cc8bc54SJed Brown *nel = 0; 4839cc8bc54SJed Brown *nen = -1; 4849cc8bc54SJed Brown *e = NULL; 4852dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 4862dde6fd4SLisandro Dalcin } 487f951a8fcSStefano Zampini 488f951a8fcSStefano Zampini /*@ 489*dce8aebaSBarry Smith DMDARestoreSubdomainCornersIS - Restores the `IS` obtained with `DMDAGetSubdomainCornersIS()` 490f951a8fcSStefano Zampini 491f951a8fcSStefano Zampini Not Collective 492f951a8fcSStefano Zampini 493d8d19677SJose E. Roman Input Parameters: 494*dce8aebaSBarry Smith + dm - the `DM` object 495f951a8fcSStefano Zampini - is - the index set 496f951a8fcSStefano Zampini 497f951a8fcSStefano Zampini Level: intermediate 498f951a8fcSStefano Zampini 499*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetSubdomainCornersIS()` 500f951a8fcSStefano Zampini @*/ 501d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm, IS *is) 502d71ae5a4SJacob Faibussowitsch { 503f951a8fcSStefano Zampini PetscFunctionBegin; 504a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 505f951a8fcSStefano Zampini PetscValidHeaderSpecific(*is, IS_CLASSID, 2); 506f951a8fcSStefano Zampini *is = NULL; 507f951a8fcSStefano Zampini PetscFunctionReturn(0); 508f951a8fcSStefano Zampini } 509