xref: /petsc/src/dm/impls/da/dagetelem.c (revision cfbed8a3daa685ef2133d2c62dbddaab837fc99c)
1454e267fSLisandro Dalcin 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>     /*I  "petscdmda.h"   I*/
3454e267fSLisandro Dalcin 
42dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_1D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
5454e267fSLisandro Dalcin {
6454e267fSLisandro Dalcin   PetscErrorCode ierr;
7454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
8454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
9454e267fSLisandro Dalcin   PetscInt       cnt=0;
105fd66863SKarl Rupp 
11454e267fSLisandro Dalcin   PetscFunctionBegin;
12454e267fSLisandro Dalcin   if (!da->e) {
13f951a8fcSStefano Zampini     PetscInt corners[2];
14f951a8fcSStefano Zampini 
1598b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
16454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,0,0,&xe,0,0);CHKERRQ(ierr);
17454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,0,0,&Xe,0,0);CHKERRQ(ierr);
18454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
19454e267fSLisandro Dalcin     da->ne = 1*(xe - xs - 1);
20854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + 2*da->ne,&da->e);CHKERRQ(ierr);
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     }
25*cfbed8a3SStefano Zampini     da->nen = 2;
26*cfbed8a3SStefano Zampini 
27f951a8fcSStefano Zampini     corners[0] = (xs  -Xs);
28f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs);
29f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,2,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
30454e267fSLisandro Dalcin   }
31454e267fSLisandro Dalcin   *nel = da->ne;
32*cfbed8a3SStefano Zampini   *nen = da->nen;
33454e267fSLisandro Dalcin   *e   = da->e;
34454e267fSLisandro Dalcin   PetscFunctionReturn(0);
35454e267fSLisandro Dalcin }
36454e267fSLisandro Dalcin 
372dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
38454e267fSLisandro Dalcin {
39454e267fSLisandro Dalcin   PetscErrorCode ierr;
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;
43*cfbed8a3SStefano Zampini   PetscInt       cnt=0, cell[4], ns=2;
44454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,
45454e267fSLisandro Dalcin                                2,3,1};
465fd66863SKarl Rupp 
47454e267fSLisandro Dalcin   PetscFunctionBegin;
48454e267fSLisandro Dalcin   if (!da->e) {
49*cfbed8a3SStefano Zampini     PetscInt corners[4],nn = 0;
50f951a8fcSStefano Zampini 
5198b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
52*cfbed8a3SStefano Zampini 
53*cfbed8a3SStefano Zampini     switch (da->elementtype) {
54*cfbed8a3SStefano Zampini     case DMDA_ELEMENT_Q1:
55*cfbed8a3SStefano Zampini       da->nen = 4;
56*cfbed8a3SStefano Zampini       break;
57*cfbed8a3SStefano Zampini     case DMDA_ELEMENT_P1:
58*cfbed8a3SStefano Zampini       da->nen = 3;
59*cfbed8a3SStefano Zampini       break;
60*cfbed8a3SStefano Zampini     default:
61*cfbed8a3SStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype);
62*cfbed8a3SStefano Zampini       break;
63*cfbed8a3SStefano Zampini     }
64*cfbed8a3SStefano Zampini     nn = da->nen;
65*cfbed8a3SStefano Zampini 
66b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2;}
67b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
68454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr);
69454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr);
70454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
71454e267fSLisandro Dalcin     ye    += ys; Ye += Ys; if (ys != Ys) ys -= 1;
72454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1);
73854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr);
74454e267fSLisandro Dalcin     for (j=ys; j<ye-1; j++) {
75454e267fSLisandro Dalcin       for (i=xs; i<xe-1; i++) {
76454e267fSLisandro Dalcin         cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs);
77454e267fSLisandro Dalcin         cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs);
78454e267fSLisandro Dalcin         cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs);
79454e267fSLisandro Dalcin         cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs);
80454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
818865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]];
82454e267fSLisandro Dalcin         }
83454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
848865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c];
85454e267fSLisandro Dalcin         }
86454e267fSLisandro Dalcin       }
87454e267fSLisandro Dalcin     }
88*cfbed8a3SStefano Zampini 
89f951a8fcSStefano Zampini     corners[0] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs);
90f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs);
91f951a8fcSStefano Zampini     corners[2] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs);
92f951a8fcSStefano Zampini     corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs);
93f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,4,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
94454e267fSLisandro Dalcin   }
95454e267fSLisandro Dalcin   *nel = da->ne;
96*cfbed8a3SStefano Zampini   *nen = da->nen;
97454e267fSLisandro Dalcin   *e   = da->e;
98454e267fSLisandro Dalcin   PetscFunctionReturn(0);
99454e267fSLisandro Dalcin }
100454e267fSLisandro Dalcin 
1012dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
102454e267fSLisandro Dalcin {
103454e267fSLisandro Dalcin   PetscErrorCode ierr;
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;
108*cfbed8a3SStefano Zampini   PetscInt       cnt=0, cell[8], ns=6;
109454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
110454e267fSLisandro Dalcin                                0,1,7,4,
111454e267fSLisandro Dalcin                                1,2,3,7,
112454e267fSLisandro Dalcin                                1,2,7,6,
113454e267fSLisandro Dalcin                                1,4,5,7,
114454e267fSLisandro Dalcin                                1,5,6,7};
1155fd66863SKarl Rupp 
116454e267fSLisandro Dalcin   PetscFunctionBegin;
117454e267fSLisandro Dalcin   if (!da->e) {
118*cfbed8a3SStefano Zampini     PetscInt corners[8],nn = 0;
119f951a8fcSStefano Zampini 
12098b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
121*cfbed8a3SStefano Zampini 
122*cfbed8a3SStefano Zampini     switch (da->elementtype) {
123*cfbed8a3SStefano Zampini     case DMDA_ELEMENT_Q1:
124*cfbed8a3SStefano Zampini       da->nen = 8;
125*cfbed8a3SStefano Zampini       break;
126*cfbed8a3SStefano Zampini     case DMDA_ELEMENT_P1:
127*cfbed8a3SStefano Zampini       da->nen = 4;
128*cfbed8a3SStefano Zampini       break;
129*cfbed8a3SStefano Zampini     default:
130*cfbed8a3SStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype);
131*cfbed8a3SStefano Zampini       break;
132*cfbed8a3SStefano Zampini     }
133*cfbed8a3SStefano Zampini     nn = da->nen;
134*cfbed8a3SStefano Zampini 
135b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6;}
136b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
137454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
138454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
139454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
140454e267fSLisandro Dalcin     ye    += ys; Ye += Ys; if (ys != Ys) ys -= 1;
141454e267fSLisandro Dalcin     ze    += zs; Ze += Zs; if (zs != Zs) zs -= 1;
142454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
143854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr);
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);
148454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
149454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
150454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
151454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
152454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
153454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
154454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(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     }
164*cfbed8a3SStefano 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);
173f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,8,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
174454e267fSLisandro Dalcin   }
175454e267fSLisandro Dalcin   *nel = da->ne;
176*cfbed8a3SStefano Zampini   *nen = da->nen;
177454e267fSLisandro Dalcin   *e   = da->e;
178454e267fSLisandro Dalcin   PetscFunctionReturn(0);
179454e267fSLisandro Dalcin }
180454e267fSLisandro Dalcin 
181f5f57ec0SBarry Smith /*@
182d4a6ed37SStefano Zampini    DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left
183d4a6ed37SStefano Zampini    corner of the non-overlapping decomposition identified by DMDAGetElements()
184d4a6ed37SStefano Zampini 
185d4a6ed37SStefano Zampini     Not Collective
186d4a6ed37SStefano Zampini 
187d4a6ed37SStefano Zampini    Input Parameter:
188d4a6ed37SStefano Zampini .     da - the DM object
189d4a6ed37SStefano Zampini 
190d4a6ed37SStefano Zampini    Output Parameters:
191d4a6ed37SStefano Zampini +     gx - the x index
192d4a6ed37SStefano Zampini .     gy - the y index
193d4a6ed37SStefano Zampini -     gz - the z index
194d4a6ed37SStefano Zampini 
195d4a6ed37SStefano Zampini    Level: intermediate
196d4a6ed37SStefano Zampini 
197d4a6ed37SStefano Zampini    Notes:
198d4a6ed37SStefano Zampini 
199d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
200d4a6ed37SStefano Zampini @*/
201d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz)
202d4a6ed37SStefano Zampini {
203d4a6ed37SStefano Zampini   PetscInt       xs,Xs;
204d4a6ed37SStefano Zampini   PetscInt       ys,Ys;
205d4a6ed37SStefano Zampini   PetscInt       zs,Zs;
206d4a6ed37SStefano Zampini   PetscBool      isda;
207d4a6ed37SStefano Zampini   PetscErrorCode ierr;
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);
214d4a6ed37SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
215d4a6ed37SStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
216d4a6ed37SStefano Zampini   ierr = DMDAGetCorners(da,&xs,&ys,&zs,NULL,NULL,NULL);CHKERRQ(ierr);
217d4a6ed37SStefano Zampini   ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr);
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 
245d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements
246d4a6ed37SStefano Zampini @*/
247d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz)
248d4a6ed37SStefano Zampini {
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   PetscErrorCode ierr;
255d4a6ed37SStefano Zampini 
256d4a6ed37SStefano Zampini   PetscFunctionBegin;
257a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
258d4a6ed37SStefano Zampini   if (mx) PetscValidIntPointer(mx,2);
259d4a6ed37SStefano Zampini   if (my) PetscValidIntPointer(my,3);
260d4a6ed37SStefano Zampini   if (mz) PetscValidIntPointer(mz,4);
261d4a6ed37SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
262d4a6ed37SStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
263d4a6ed37SStefano Zampini   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
264d4a6ed37SStefano Zampini   ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr);
265d4a6ed37SStefano Zampini   xe  += xs; if (xs != Xs) xs -= 1;
266d4a6ed37SStefano Zampini   ye  += ys; if (ys != Ys) ys -= 1;
267d4a6ed37SStefano Zampini   ze  += zs; if (zs != Zs) zs -= 1;
268d4a6ed37SStefano Zampini   if (mx) *mx  = 0;
269d4a6ed37SStefano Zampini   if (my) *my  = 0;
270d4a6ed37SStefano Zampini   if (mz) *mz  = 0;
271d4a6ed37SStefano Zampini   ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
272d4a6ed37SStefano Zampini   switch (dim) {
273d4a6ed37SStefano Zampini   case 3:
274d4a6ed37SStefano Zampini     if (mz) *mz = ze - zs - 1;
275d4a6ed37SStefano Zampini   case 2:
276d4a6ed37SStefano Zampini     if (my) *my = ye - ys - 1;
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 /*@
2852dde6fd4SLisandro Dalcin       DMDASetElementType - Sets the element type to be returned by DMDAGetElements()
2862dde6fd4SLisandro Dalcin 
2872dde6fd4SLisandro Dalcin     Not Collective
2882dde6fd4SLisandro Dalcin 
2892dde6fd4SLisandro Dalcin    Input Parameter:
2902dde6fd4SLisandro Dalcin .     da - the DMDA object
2912dde6fd4SLisandro Dalcin 
2922dde6fd4SLisandro Dalcin    Output Parameters:
2936420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
2942dde6fd4SLisandro Dalcin 
2952dde6fd4SLisandro Dalcin    Level: intermediate
2962dde6fd4SLisandro Dalcin 
2972dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
2982dde6fd4SLisandro Dalcin @*/
2992dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
3002dde6fd4SLisandro Dalcin {
3012dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
3022dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
303f951a8fcSStefano Zampini   PetscBool      isda;
3042dde6fd4SLisandro Dalcin 
3052dde6fd4SLisandro Dalcin   PetscFunctionBegin;
306a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
3072dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
308f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
309f951a8fcSStefano Zampini   if (!isda) PetscFunctionReturn(0);
3102dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
3112dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
312f951a8fcSStefano Zampini     ierr = ISDestroy(&dd->ecorners);CHKERRQ(ierr);
3138865f1eaSKarl Rupp 
3142dde6fd4SLisandro Dalcin     dd->elementtype = etype;
3152dde6fd4SLisandro Dalcin     dd->ne          = 0;
316*cfbed8a3SStefano Zampini     dd->nen         = 0;
3170298fd71SBarry Smith     dd->e           = NULL;
3182dde6fd4SLisandro Dalcin   }
3192dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
3202dde6fd4SLisandro Dalcin }
3212dde6fd4SLisandro Dalcin 
322f5f57ec0SBarry Smith /*@
3232dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
3242dde6fd4SLisandro Dalcin 
3252dde6fd4SLisandro Dalcin     Not Collective
3262dde6fd4SLisandro Dalcin 
3272dde6fd4SLisandro Dalcin    Input Parameter:
3282dde6fd4SLisandro Dalcin .     da - the DMDA object
3292dde6fd4SLisandro Dalcin 
3302dde6fd4SLisandro Dalcin    Output Parameters:
3316420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
3322dde6fd4SLisandro Dalcin 
3332dde6fd4SLisandro Dalcin    Level: intermediate
3342dde6fd4SLisandro Dalcin 
3352dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
3362dde6fd4SLisandro Dalcin @*/
3372dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
3382dde6fd4SLisandro Dalcin {
3392dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
340f951a8fcSStefano Zampini   PetscErrorCode ierr;
341f951a8fcSStefano Zampini   PetscBool      isda;
3425fd66863SKarl Rupp 
3432dde6fd4SLisandro Dalcin   PetscFunctionBegin;
344a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
3452dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
346f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
347f951a8fcSStefano Zampini   if (!isda) SETERRQ1(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 
37793386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
3782dde6fd4SLisandro Dalcin @*/
3792dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
380454e267fSLisandro Dalcin {
381c73cfb54SMatthew G. Knepley   PetscInt       dim;
382454e267fSLisandro Dalcin   PetscErrorCode ierr;
38399c57625SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
384f951a8fcSStefano Zampini   PetscBool      isda;
3855fd66863SKarl Rupp 
386454e267fSLisandro Dalcin   PetscFunctionBegin;
387a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA);
388f951a8fcSStefano Zampini   PetscValidIntPointer(nel,2);
389f951a8fcSStefano Zampini   PetscValidIntPointer(nen,3);
390f951a8fcSStefano Zampini   PetscValidPointer(e,4);
391f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
392f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
393f951a8fcSStefano Zampini   if (dd->stencil_type == DMDA_STENCIL_STAR) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMDAGetElements() requires you use a stencil type of DMDA_STENCIL_BOX");
394c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
395*cfbed8a3SStefano Zampini   if (dd->e) {
396*cfbed8a3SStefano Zampini     *nel = dd->ne;
397*cfbed8a3SStefano Zampini     *nen = dd->nen;
398*cfbed8a3SStefano Zampini     *e   = dd->e;
399*cfbed8a3SStefano Zampini     PetscFunctionReturn(0);
400*cfbed8a3SStefano Zampini   }
401c73cfb54SMatthew G. Knepley   if (dim==-1) {
4020298fd71SBarry Smith     *nel = 0; *nen = 0; *e = NULL;
403c73cfb54SMatthew G. Knepley   } else if (dim==1) {
4042dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
405c73cfb54SMatthew G. Knepley   } else if (dim==2) {
4062dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
407c73cfb54SMatthew G. Knepley   } else if (dim==3) {
4082dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
409c73cfb54SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",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)
415f951a8fcSStefano Zampini                                  of the non-overlapping decomposition identified by DMDAGetElements
416f951a8fcSStefano Zampini 
417f951a8fcSStefano Zampini     Not Collective
418f951a8fcSStefano Zampini 
419f951a8fcSStefano Zampini    Input Parameter:
420f951a8fcSStefano Zampini .     dm - the DM object
421f951a8fcSStefano Zampini 
422f951a8fcSStefano Zampini    Output Parameters:
423f951a8fcSStefano Zampini .     is - the index set
424f951a8fcSStefano Zampini 
425f951a8fcSStefano Zampini    Level: intermediate
426f951a8fcSStefano Zampini 
42795452b02SPatrick Sanan    Notes:
42895452b02SPatrick Sanan     Call DMDARestoreSubdomainCornersIS() once you have finished accessing the index set.
429f951a8fcSStefano Zampini 
430f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElementsCornersIS()
431f951a8fcSStefano Zampini @*/
432d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetSubdomainCornersIS(DM dm,IS *is)
433f951a8fcSStefano Zampini {
434f951a8fcSStefano Zampini   PetscErrorCode ierr;
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);
441f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
442f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
443f951a8fcSStefano Zampini   if (dd->stencil_type == DMDA_STENCIL_STAR) SETERRQ(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 
448f951a8fcSStefano Zampini     ierr = DMDAGetElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
449f951a8fcSStefano Zampini     ierr = DMDARestoreElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
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 
4602dde6fd4SLisandro Dalcin    Input Parameter:
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 
4742dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
4752dde6fd4SLisandro Dalcin @*/
4762dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
4772dde6fd4SLisandro Dalcin {
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 
494f951a8fcSStefano Zampini    Input Parameter:
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 
502d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetSubdomainCornersIS()
503f951a8fcSStefano Zampini @*/
504d4a6ed37SStefano Zampini PetscErrorCode  DMDARestoreSubdomainCornersIS(DM dm,IS *is)
505f951a8fcSStefano Zampini {
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