xref: /petsc/src/dm/impls/da/dagetelem.c (revision 6d8694c4fbab79f9439f1ad13c0386ba7ee1ca4b)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I  "petscdmda.h"   I*/
2454e267fSLisandro Dalcin 
DMDAGetElements_1D(DM dm,PetscInt * nel,PetscInt * nen,const PetscInt * e[])3d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_1D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
4d71ae5a4SJacob Faibussowitsch {
5454e267fSLisandro Dalcin   DM_DA   *da = (DM_DA *)dm->data;
6454e267fSLisandro Dalcin   PetscInt i, xs, xe, Xs, Xe;
7454e267fSLisandro Dalcin   PetscInt cnt = 0;
85fd66863SKarl Rupp 
9454e267fSLisandro Dalcin   PetscFunctionBegin;
10454e267fSLisandro Dalcin   if (!da->e) {
11f951a8fcSStefano Zampini     PetscInt corners[2];
12f951a8fcSStefano Zampini 
137a8be351SBarry Smith     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
149566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xe, NULL, NULL));
159566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(dm, &Xs, NULL, NULL, &Xe, NULL, NULL));
169371c9d4SSatish Balay     xe += xs;
179371c9d4SSatish Balay     Xe += Xs;
189371c9d4SSatish Balay     if (xs != Xs) xs -= 1;
19454e267fSLisandro Dalcin     da->ne = 1 * (xe - xs - 1);
209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + 2 * da->ne, &da->e));
21454e267fSLisandro Dalcin     for (i = xs; i < xe - 1; i++) {
22454e267fSLisandro Dalcin       da->e[cnt++] = (i - Xs);
23454e267fSLisandro Dalcin       da->e[cnt++] = (i - Xs + 1);
24454e267fSLisandro Dalcin     }
25cfbed8a3SStefano Zampini     da->nen = 2;
26cfbed8a3SStefano Zampini 
27f951a8fcSStefano Zampini     corners[0] = (xs - Xs);
28f951a8fcSStefano Zampini     corners[1] = (xe - 1 - Xs);
299566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, corners, PETSC_COPY_VALUES, &da->ecorners));
30454e267fSLisandro Dalcin   }
31454e267fSLisandro Dalcin   *nel = da->ne;
32cfbed8a3SStefano Zampini   *nen = da->nen;
33454e267fSLisandro Dalcin   *e   = da->e;
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35454e267fSLisandro Dalcin }
36454e267fSLisandro Dalcin 
DMDAGetElements_2D(DM dm,PetscInt * nel,PetscInt * nen,const PetscInt * e[])37d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_2D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
38d71ae5a4SJacob Faibussowitsch {
39454e267fSLisandro Dalcin   DM_DA   *da = (DM_DA *)dm->data;
40454e267fSLisandro Dalcin   PetscInt i, xs, xe, Xs, Xe;
41454e267fSLisandro Dalcin   PetscInt j, ys, ye, Ys, Ye;
42cfbed8a3SStefano Zampini   PetscInt cnt = 0, cell[4], ns = 2;
439371c9d4SSatish Balay   PetscInt c, split[] = {0, 1, 3, 2, 3, 1};
445fd66863SKarl Rupp 
45454e267fSLisandro Dalcin   PetscFunctionBegin;
46454e267fSLisandro Dalcin   if (!da->e) {
47cfbed8a3SStefano Zampini     PetscInt corners[4], nn = 0;
48f951a8fcSStefano Zampini 
497a8be351SBarry Smith     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
50cfbed8a3SStefano Zampini 
51cfbed8a3SStefano Zampini     switch (da->elementtype) {
52d71ae5a4SJacob Faibussowitsch     case DMDA_ELEMENT_Q1:
53d71ae5a4SJacob Faibussowitsch       da->nen = 4;
54d71ae5a4SJacob Faibussowitsch       break;
55d71ae5a4SJacob Faibussowitsch     case DMDA_ELEMENT_P1:
56d71ae5a4SJacob Faibussowitsch       da->nen = 3;
57d71ae5a4SJacob Faibussowitsch       break;
58d71ae5a4SJacob Faibussowitsch     default:
59d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype);
60cfbed8a3SStefano Zampini     }
61cfbed8a3SStefano Zampini     nn = da->nen;
62cfbed8a3SStefano Zampini 
63ad540459SPierre Jolivet     if (da->elementtype == DMDA_ELEMENT_P1) ns = 2;
64ad540459SPierre Jolivet     if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1;
659566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, &ys, NULL, &xe, &ye, NULL));
669566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, NULL, &Xe, &Ye, NULL));
679371c9d4SSatish Balay     xe += xs;
689371c9d4SSatish Balay     Xe += Xs;
699371c9d4SSatish Balay     if (xs != Xs) xs -= 1;
709371c9d4SSatish Balay     ye += ys;
719371c9d4SSatish Balay     Ye += Ys;
729371c9d4SSatish Balay     if (ys != Ys) ys -= 1;
73454e267fSLisandro Dalcin     da->ne = ns * (xe - xs - 1) * (ye - ys - 1);
749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e));
75454e267fSLisandro Dalcin     for (j = ys; j < ye - 1; j++) {
76454e267fSLisandro Dalcin       for (i = xs; i < xe - 1; i++) {
77454e267fSLisandro Dalcin         cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs);
78454e267fSLisandro Dalcin         cell[1] = (i - Xs + 1) + (j - Ys) * (Xe - Xs);
79454e267fSLisandro Dalcin         cell[2] = (i - Xs + 1) + (j - Ys + 1) * (Xe - Xs);
80454e267fSLisandro Dalcin         cell[3] = (i - Xs) + (j - Ys + 1) * (Xe - Xs);
81454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
828865f1eaSKarl Rupp           for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]];
83454e267fSLisandro Dalcin         }
84454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
858865f1eaSKarl Rupp           for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c];
86454e267fSLisandro Dalcin         }
87454e267fSLisandro Dalcin       }
88454e267fSLisandro Dalcin     }
89cfbed8a3SStefano Zampini 
90f951a8fcSStefano Zampini     corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs);
91f951a8fcSStefano Zampini     corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs);
92f951a8fcSStefano Zampini     corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs);
93f951a8fcSStefano Zampini     corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs);
949566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 4, corners, PETSC_COPY_VALUES, &da->ecorners));
95454e267fSLisandro Dalcin   }
96454e267fSLisandro Dalcin   *nel = da->ne;
97cfbed8a3SStefano Zampini   *nen = da->nen;
98454e267fSLisandro Dalcin   *e   = da->e;
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100454e267fSLisandro Dalcin }
101454e267fSLisandro Dalcin 
DMDAGetElements_3D(DM dm,PetscInt * nel,PetscInt * nen,const PetscInt * e[])102d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetElements_3D(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
103d71ae5a4SJacob Faibussowitsch {
104454e267fSLisandro Dalcin   DM_DA   *da = (DM_DA *)dm->data;
105454e267fSLisandro Dalcin   PetscInt i, xs, xe, Xs, Xe;
106454e267fSLisandro Dalcin   PetscInt j, ys, ye, Ys, Ye;
107454e267fSLisandro Dalcin   PetscInt k, zs, ze, Zs, Ze;
108cfbed8a3SStefano Zampini   PetscInt cnt = 0, cell[8], ns = 6;
1099371c9d4SSatish 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};
1105fd66863SKarl Rupp 
111454e267fSLisandro Dalcin   PetscFunctionBegin;
112454e267fSLisandro Dalcin   if (!da->e) {
113cfbed8a3SStefano Zampini     PetscInt corners[8], nn = 0;
114f951a8fcSStefano Zampini 
1157a8be351SBarry Smith     PetscCheck(da->s, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot get elements for DMDA with zero stencil width");
116cfbed8a3SStefano Zampini 
117cfbed8a3SStefano Zampini     switch (da->elementtype) {
118d71ae5a4SJacob Faibussowitsch     case DMDA_ELEMENT_Q1:
119d71ae5a4SJacob Faibussowitsch       da->nen = 8;
120d71ae5a4SJacob Faibussowitsch       break;
121d71ae5a4SJacob Faibussowitsch     case DMDA_ELEMENT_P1:
122d71ae5a4SJacob Faibussowitsch       da->nen = 4;
123d71ae5a4SJacob Faibussowitsch       break;
124d71ae5a4SJacob Faibussowitsch     default:
125d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown element type %d", da->elementtype);
126cfbed8a3SStefano Zampini     }
127cfbed8a3SStefano Zampini     nn = da->nen;
128cfbed8a3SStefano Zampini 
129ad540459SPierre Jolivet     if (da->elementtype == DMDA_ELEMENT_P1) ns = 6;
130ad540459SPierre Jolivet     if (da->elementtype == DMDA_ELEMENT_Q1) ns = 1;
1319566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, &ys, &zs, &xe, &ye, &ze));
1329566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(dm, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze));
1339371c9d4SSatish Balay     xe += xs;
1349371c9d4SSatish Balay     Xe += Xs;
1359371c9d4SSatish Balay     if (xs != Xs) xs -= 1;
1369371c9d4SSatish Balay     ye += ys;
1379371c9d4SSatish Balay     Ye += Ys;
1389371c9d4SSatish Balay     if (ys != Ys) ys -= 1;
1399371c9d4SSatish Balay     ze += zs;
1409371c9d4SSatish Balay     Ze += Zs;
1419371c9d4SSatish Balay     if (zs != Zs) zs -= 1;
142454e267fSLisandro Dalcin     da->ne = ns * (xe - xs - 1) * (ye - ys - 1) * (ze - zs - 1);
1439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + nn * da->ne, &da->e));
144454e267fSLisandro Dalcin     for (k = zs; k < ze - 1; k++) {
145454e267fSLisandro Dalcin       for (j = ys; j < ye - 1; j++) {
146454e267fSLisandro Dalcin         for (i = xs; i < xe - 1; i++) {
147454e267fSLisandro Dalcin           cell[0] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
1482da392ccSBarry Smith           cell[1] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
1492da392ccSBarry Smith           cell[2] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
1502da392ccSBarry Smith           cell[3] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k - Zs) * (Xe - Xs) * (Ye - Ys);
1512da392ccSBarry Smith           cell[4] = (i - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
1522da392ccSBarry Smith           cell[5] = (i + 1 - Xs) + (j - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
1532da392ccSBarry Smith           cell[6] = (i + 1 - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
1542da392ccSBarry Smith           cell[7] = (i - Xs) + (j + 1 - Ys) * (Xe - Xs) + (k + 1 - Zs) * (Xe - Xs) * (Ye - Ys);
155454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
1568865f1eaSKarl Rupp             for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[split[c]];
157454e267fSLisandro Dalcin           }
158454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
1598865f1eaSKarl Rupp             for (c = 0; c < ns * nn; c++) da->e[cnt++] = cell[c];
160454e267fSLisandro Dalcin           }
161454e267fSLisandro Dalcin         }
162454e267fSLisandro Dalcin       }
163454e267fSLisandro Dalcin     }
164cfbed8a3SStefano Zampini 
165f951a8fcSStefano Zampini     corners[0] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
166f951a8fcSStefano Zampini     corners[1] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
167f951a8fcSStefano Zampini     corners[2] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
168f951a8fcSStefano Zampini     corners[3] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
169f951a8fcSStefano Zampini     corners[4] = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
170f951a8fcSStefano Zampini     corners[5] = (xe - 1 - Xs) + (ys - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
171f951a8fcSStefano Zampini     corners[6] = (xs - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
172f951a8fcSStefano Zampini     corners[7] = (xe - 1 - Xs) + (ye - 1 - Ys) * (Xe - Xs) + (ze - 1 - Zs) * (Xe - Xs) * (Ye - Ys);
1739566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 8, corners, PETSC_COPY_VALUES, &da->ecorners));
174454e267fSLisandro Dalcin   }
175454e267fSLisandro Dalcin   *nel = da->ne;
176cfbed8a3SStefano Zampini   *nen = da->nen;
177454e267fSLisandro Dalcin   *e   = da->e;
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179454e267fSLisandro Dalcin }
180454e267fSLisandro Dalcin 
181f5f57ec0SBarry Smith /*@
18272fd0fbdSBarry Smith   DMDAGetElementsCorners - Returns the global (i,j,k) indices of the lower left
18372fd0fbdSBarry Smith   corner of the non-overlapping decomposition of elements identified by `DMDAGetElements()`
184d4a6ed37SStefano Zampini 
185d4a6ed37SStefano Zampini   Not Collective
186d4a6ed37SStefano Zampini 
187d4a6ed37SStefano Zampini   Input Parameter:
188dce8aebaSBarry Smith . da - the `DMDA` object
189d4a6ed37SStefano Zampini 
190d4a6ed37SStefano Zampini   Output Parameters:
19172fd0fbdSBarry Smith + gx - the i index
19272fd0fbdSBarry Smith . gy - the j index
19372fd0fbdSBarry Smith - gz - the k index
194d4a6ed37SStefano Zampini 
195d4a6ed37SStefano Zampini   Level: intermediate
196d4a6ed37SStefano Zampini 
19712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDAGetElementsSizes()`,
19812b4a537SBarry Smith           `DMDAGetElementsCornersIS()`, `DMDARestoreElementsCornersIS()`
199d4a6ed37SStefano Zampini @*/
DMDAGetElementsCorners(DM da,PeOp PetscInt * gx,PeOp PetscInt * gy,PeOp PetscInt * gz)200*ce78bad3SBarry Smith PetscErrorCode DMDAGetElementsCorners(DM da, PeOp PetscInt *gx, PeOp PetscInt *gy, PeOp 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);
2094f572ea9SToby Isaac   if (gx) PetscAssertPointer(gx, 2);
2104f572ea9SToby Isaac   if (gy) PetscAssertPointer(gy, 3);
2114f572ea9SToby Isaac   if (gz) PetscAssertPointer(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;
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
223d4a6ed37SStefano Zampini }
224d4a6ed37SStefano Zampini 
225d4a6ed37SStefano Zampini /*@
22612b4a537SBarry Smith   DMDAGetElementsSizes - Gets the local number of elements per coordinate direction for the non-overlapping decomposition identified by `DMDAGetElements()`
227d4a6ed37SStefano Zampini 
228d4a6ed37SStefano Zampini   Not Collective
229d4a6ed37SStefano Zampini 
230d4a6ed37SStefano Zampini   Input Parameter:
231dce8aebaSBarry 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 
240dce8aebaSBarry Smith   Note:
24172fd0fbdSBarry Smith   Returns the same number of elements, irrespective of the `DMDAElementType`
242d4a6ed37SStefano Zampini 
24312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDAGetElementsCorners()`
244d4a6ed37SStefano Zampini @*/
DMDAGetElementsSizes(DM da,PeOp PetscInt * mx,PeOp PetscInt * my,PeOp PetscInt * mz)245*ce78bad3SBarry Smith PetscErrorCode DMDAGetElementsSizes(DM da, PeOp PetscInt *mx, PeOp PetscInt *my, PeOp 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);
2554f572ea9SToby Isaac   if (mx) PetscAssertPointer(mx, 2);
2564f572ea9SToby Isaac   if (my) PetscAssertPointer(my, 3);
2574f572ea9SToby Isaac   if (mz) PetscAssertPointer(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   }
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
282d4a6ed37SStefano Zampini }
283d4a6ed37SStefano Zampini 
284d4a6ed37SStefano Zampini /*@
285dce8aebaSBarry Smith   DMDASetElementType - Sets the element type to be returned by `DMDAGetElements()`
2862dde6fd4SLisandro Dalcin 
2872dde6fd4SLisandro Dalcin   Not Collective
2882dde6fd4SLisandro Dalcin 
2892dde6fd4SLisandro Dalcin   Input Parameter:
290dce8aebaSBarry Smith . da - the `DMDA` object
2912dde6fd4SLisandro Dalcin 
2922fe279fdSBarry Smith   Output Parameter:
293dce8aebaSBarry Smith . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1`
2942dde6fd4SLisandro Dalcin 
2952dde6fd4SLisandro Dalcin   Level: intermediate
2962dde6fd4SLisandro Dalcin 
29712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDAGetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()`,
29812b4a537SBarry Smith           `DMDA_ELEMENT_P1`, `DMDA_ELEMENT_Q1`
2992dde6fd4SLisandro Dalcin @*/
DMDASetElementType(DM da,DMDAElementType etype)300d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype)
301d71ae5a4SJacob Faibussowitsch {
3022dde6fd4SLisandro Dalcin   DM_DA    *dd = (DM_DA *)da->data;
303f951a8fcSStefano Zampini   PetscBool isda;
3042dde6fd4SLisandro Dalcin 
3052dde6fd4SLisandro Dalcin   PetscFunctionBegin;
306a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
3072dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, etype, 2);
3089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)da, DMDA, &isda));
3093ba16761SJacob Faibussowitsch   if (!isda) PetscFunctionReturn(PETSC_SUCCESS);
3102dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
3119566063dSJacob Faibussowitsch     PetscCall(PetscFree(dd->e));
3129566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&dd->ecorners));
3138865f1eaSKarl Rupp 
3142dde6fd4SLisandro Dalcin     dd->elementtype = etype;
3152dde6fd4SLisandro Dalcin     dd->ne          = 0;
316cfbed8a3SStefano Zampini     dd->nen         = 0;
3170298fd71SBarry Smith     dd->e           = NULL;
3182dde6fd4SLisandro Dalcin   }
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3202dde6fd4SLisandro Dalcin }
3212dde6fd4SLisandro Dalcin 
322f5f57ec0SBarry Smith /*@
323dce8aebaSBarry Smith   DMDAGetElementType - Gets the element type to be returned by `DMDAGetElements()`
3242dde6fd4SLisandro Dalcin 
3252dde6fd4SLisandro Dalcin   Not Collective
3262dde6fd4SLisandro Dalcin 
3272dde6fd4SLisandro Dalcin   Input Parameter:
328dce8aebaSBarry Smith . da - the `DMDA` object
3292dde6fd4SLisandro Dalcin 
3302fe279fdSBarry Smith   Output Parameter:
331dce8aebaSBarry Smith . etype - the element type, currently either `DMDA_ELEMENT_P1` or `DMDA_ELEMENT_Q1`
3322dde6fd4SLisandro Dalcin 
3332dde6fd4SLisandro Dalcin   Level: intermediate
3342dde6fd4SLisandro Dalcin 
33512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElements()`,
33612b4a537SBarry Smith           `DMDA_ELEMENT_P1`, `DMDA_ELEMENT_Q1`
3372dde6fd4SLisandro Dalcin @*/
DMDAGetElementType(DM da,DMDAElementType * etype)338d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype)
339d71ae5a4SJacob 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);
3454f572ea9SToby Isaac   PetscAssertPointer(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;
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3502dde6fd4SLisandro Dalcin }
3512dde6fd4SLisandro Dalcin 
3522dde6fd4SLisandro Dalcin /*@C
35372fd0fbdSBarry Smith   DMDAGetElements - Gets an array containing the indices (in local indexing)
3542dde6fd4SLisandro Dalcin   of all the local elements
3552dde6fd4SLisandro Dalcin 
356659f25fdSBarry Smith   Not Collective
3572dde6fd4SLisandro Dalcin 
3582dde6fd4SLisandro Dalcin   Input Parameter:
359dce8aebaSBarry Smith . dm - the `DMDA` object
3602dde6fd4SLisandro Dalcin 
3612dde6fd4SLisandro Dalcin   Output Parameters:
3622dde6fd4SLisandro Dalcin + nel - number of local elements
36312b4a537SBarry Smith . nen - number of nodes in each element (for example in one dimension it is 2, in two dimensions it is 3 (for `DMDA_ELEMENT_P1`) and 4
36412b4a537SBarry Smith         (for `DMDA_ELEMENT_Q1`)
365f13dfd9eSBarry Smith - e   - the local indices of the elements' vertices, of length `nel` * `nen`
3662dde6fd4SLisandro Dalcin 
3672dde6fd4SLisandro Dalcin   Level: intermediate
3682dde6fd4SLisandro Dalcin 
36993386343SJed Brown   Notes:
370dce8aebaSBarry Smith   Call `DMDARestoreElements()` once you have finished accessing the elements.
37193386343SJed Brown 
37245950774SSatish Balay   Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.
373009f0576SBarry Smith 
374dce8aebaSBarry 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.
375009f0576SBarry Smith 
37612b4a537SBarry Smith   Fortran Note:
377659f25fdSBarry Smith   Use
378659f25fdSBarry Smith .vb
379659f25fdSBarry Smith    PetscScalar, pointer :: e(:)
380659f25fdSBarry Smith .ve
381659f25fdSBarry Smith   to declare the element array
382659f25fdSBarry Smith 
38312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`,
38412b4a537SBarry Smith           `DMGlobalToLocalBegin()`, `DMLocalToGlobalBegin()`, `DMDARestoreElements()`, `DMDA_ELEMENT_P1`, `DMDA_ELEMENT_Q1`, `DMDAGetElementsSizes()`,
38512b4a537SBarry Smith           `DMDAGetElementsCorners()`
3862dde6fd4SLisandro Dalcin @*/
DMDAGetElements(DM dm,PetscInt * nel,PetscInt * nen,const PetscInt * e[])387d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
388d71ae5a4SJacob Faibussowitsch {
389c73cfb54SMatthew G. Knepley   PetscInt  dim;
39099c57625SBarry Smith   DM_DA    *dd = (DM_DA *)dm->data;
391f951a8fcSStefano Zampini   PetscBool isda;
3925fd66863SKarl Rupp 
393454e267fSLisandro Dalcin   PetscFunctionBegin;
394a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
3954f572ea9SToby Isaac   PetscAssertPointer(nel, 2);
3964f572ea9SToby Isaac   PetscAssertPointer(nen, 3);
3974f572ea9SToby Isaac   PetscAssertPointer(e, 4);
3989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
3997a8be351SBarry Smith   PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name);
40008401ef6SPierre 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");
4019566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
402cfbed8a3SStefano Zampini   if (dd->e) {
403cfbed8a3SStefano Zampini     *nel = dd->ne;
404cfbed8a3SStefano Zampini     *nen = dd->nen;
405cfbed8a3SStefano Zampini     *e   = dd->e;
4063ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
407cfbed8a3SStefano Zampini   }
408c73cfb54SMatthew G. Knepley   if (dim == -1) {
4099371c9d4SSatish Balay     *nel = 0;
4109371c9d4SSatish Balay     *nen = 0;
4119371c9d4SSatish Balay     *e   = NULL;
412c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
4139566063dSJacob Faibussowitsch     PetscCall(DMDAGetElements_1D(dm, nel, nen, e));
414c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
4159566063dSJacob Faibussowitsch     PetscCall(DMDAGetElements_2D(dm, nel, nen, e));
416c73cfb54SMatthew G. Knepley   } else if (dim == 3) {
4179566063dSJacob Faibussowitsch     PetscCall(DMDAGetElements_3D(dm, nel, nen, e));
41863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
420454e267fSLisandro Dalcin }
4212dde6fd4SLisandro Dalcin 
422f951a8fcSStefano Zampini /*@
42372fd0fbdSBarry Smith   DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local indexing)
424dce8aebaSBarry Smith   of the non-overlapping decomposition identified by `DMDAGetElements()`
425f951a8fcSStefano Zampini 
426f951a8fcSStefano Zampini   Not Collective
427f951a8fcSStefano Zampini 
428f951a8fcSStefano Zampini   Input Parameter:
429dce8aebaSBarry Smith . dm - the `DMDA` object
430f951a8fcSStefano Zampini 
4312fe279fdSBarry Smith   Output Parameter:
432f951a8fcSStefano Zampini . is - the index set
433f951a8fcSStefano Zampini 
434f951a8fcSStefano Zampini   Level: intermediate
435f951a8fcSStefano Zampini 
436dce8aebaSBarry Smith   Note:
437dce8aebaSBarry Smith   Call `DMDARestoreSubdomainCornersIS()` once you have finished accessing the index set.
438f951a8fcSStefano Zampini 
43912b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`, `DMDARestoreElementsCornersIS()`,
44012b4a537SBarry Smith           `DMDAGetElementsSizes()`, `DMDAGetElementsCorners()`
441f951a8fcSStefano Zampini @*/
DMDAGetSubdomainCornersIS(DM dm,IS * is)442d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetSubdomainCornersIS(DM dm, IS *is)
443d71ae5a4SJacob Faibussowitsch {
444f951a8fcSStefano Zampini   DM_DA    *dd = (DM_DA *)dm->data;
445f951a8fcSStefano Zampini   PetscBool isda;
446f951a8fcSStefano Zampini 
447f951a8fcSStefano Zampini   PetscFunctionBegin;
448a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
4494f572ea9SToby Isaac   PetscAssertPointer(is, 2);
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
4517a8be351SBarry Smith   PetscCheck(isda, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for DM type %s", ((PetscObject)dm)->type_name);
45208401ef6SPierre 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");
453f951a8fcSStefano Zampini   if (!dd->ecorners) { /* compute elements if not yet done */
454f951a8fcSStefano Zampini     const PetscInt *e;
455f951a8fcSStefano Zampini     PetscInt        nel, nen;
456f951a8fcSStefano Zampini 
4579566063dSJacob Faibussowitsch     PetscCall(DMDAGetElements(dm, &nel, &nen, &e));
4589566063dSJacob Faibussowitsch     PetscCall(DMDARestoreElements(dm, &nel, &nen, &e));
459f951a8fcSStefano Zampini   }
460f951a8fcSStefano Zampini   *is = dd->ecorners;
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
462f951a8fcSStefano Zampini }
463f951a8fcSStefano Zampini 
4642dde6fd4SLisandro Dalcin /*@C
465dce8aebaSBarry Smith   DMDARestoreElements - Restores the array obtained with `DMDAGetElements()`
4662dde6fd4SLisandro Dalcin 
467659f25fdSBarry Smith   Not Collective
4682dde6fd4SLisandro Dalcin 
469d8d19677SJose E. Roman   Input Parameters:
47020f4b53cSBarry Smith + dm  - the `DM` object
4712dde6fd4SLisandro Dalcin . nel - number of local elements
472659f25fdSBarry Smith . nen - number of nodes in each element
473c2cd2aa3SJed Brown - e   - the local indices of the elements' vertices
4742dde6fd4SLisandro Dalcin 
4752dde6fd4SLisandro Dalcin   Level: intermediate
4762dde6fd4SLisandro Dalcin 
477dce8aebaSBarry Smith   Note:
478dce8aebaSBarry Smith   This restore signals the `DMDA` object that you no longer need access to the array information.
479009f0576SBarry Smith 
48012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetElements()`
4812dde6fd4SLisandro Dalcin @*/
DMDARestoreElements(DM dm,PetscInt * nel,PetscInt * nen,const PetscInt * e[])482d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreElements(DM dm, PetscInt *nel, PetscInt *nen, const PetscInt *e[])
483d71ae5a4SJacob Faibussowitsch {
4842dde6fd4SLisandro Dalcin   PetscFunctionBegin;
485a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
4864f572ea9SToby Isaac   PetscAssertPointer(nel, 2);
4874f572ea9SToby Isaac   PetscAssertPointer(nen, 3);
4884f572ea9SToby Isaac   PetscAssertPointer(e, 4);
489659f25fdSBarry Smith   if (nel) *nel = 0;
490659f25fdSBarry Smith   if (nen) *nen = -1;
491659f25fdSBarry Smith   if (e) *e = NULL;
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4932dde6fd4SLisandro Dalcin }
494f951a8fcSStefano Zampini 
495f951a8fcSStefano Zampini /*@
496dce8aebaSBarry Smith   DMDARestoreSubdomainCornersIS - Restores the `IS` obtained with `DMDAGetSubdomainCornersIS()`
497f951a8fcSStefano Zampini 
498f951a8fcSStefano Zampini   Not Collective
499f951a8fcSStefano Zampini 
500d8d19677SJose E. Roman   Input Parameters:
501dce8aebaSBarry Smith + dm - the `DM` object
502f951a8fcSStefano Zampini - is - the index set
503f951a8fcSStefano Zampini 
504f951a8fcSStefano Zampini   Level: intermediate
505f951a8fcSStefano Zampini 
50612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAElementType`, `DMDASetElementType()`, `DMDAGetSubdomainCornersIS()`
507f951a8fcSStefano Zampini @*/
DMDARestoreSubdomainCornersIS(DM dm,IS * is)508d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm, IS *is)
509d71ae5a4SJacob Faibussowitsch {
510f951a8fcSStefano Zampini   PetscFunctionBegin;
511a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
512f951a8fcSStefano Zampini   PetscValidHeaderSpecific(*is, IS_CLASSID, 2);
513f951a8fcSStefano Zampini   *is = NULL;
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515f951a8fcSStefano Zampini }
516