xref: /petsc/src/dm/impls/da/dagetelem.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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     }
25f951a8fcSStefano Zampini     corners[0] = (xs  -Xs);
26f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs);
27f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,2,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
28454e267fSLisandro Dalcin   }
29454e267fSLisandro Dalcin   *nel = da->ne;
30454e267fSLisandro Dalcin   *nen = 2;
31454e267fSLisandro Dalcin   *e   = da->e;
32454e267fSLisandro Dalcin   PetscFunctionReturn(0);
33454e267fSLisandro Dalcin }
34454e267fSLisandro Dalcin 
352dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
36454e267fSLisandro Dalcin {
37454e267fSLisandro Dalcin   PetscErrorCode ierr;
38454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
39454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
40454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
41454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[4], ns=2, nn=3;
42454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,
43454e267fSLisandro Dalcin                                2,3,1};
445fd66863SKarl Rupp 
45454e267fSLisandro Dalcin   PetscFunctionBegin;
46b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_P1) {nn=3;}
47b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_Q1) {nn=4;}
48454e267fSLisandro Dalcin   if (!da->e) {
49f951a8fcSStefano Zampini     PetscInt corners[4];
50f951a8fcSStefano Zampini 
5198b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
52b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2;}
53b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
54454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr);
55454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr);
56454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
57454e267fSLisandro Dalcin     ye    += ys; Ye += Ys; if (ys != Ys) ys -= 1;
58454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1);
59854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr);
60454e267fSLisandro Dalcin     for (j=ys; j<ye-1; j++) {
61454e267fSLisandro Dalcin       for (i=xs; i<xe-1; i++) {
62454e267fSLisandro Dalcin         cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs);
63454e267fSLisandro Dalcin         cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs);
64454e267fSLisandro Dalcin         cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs);
65454e267fSLisandro Dalcin         cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs);
66454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
678865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]];
68454e267fSLisandro Dalcin         }
69454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
708865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c];
71454e267fSLisandro Dalcin         }
72454e267fSLisandro Dalcin       }
73454e267fSLisandro Dalcin     }
74f951a8fcSStefano Zampini     corners[0] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs);
75f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs);
76f951a8fcSStefano Zampini     corners[2] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs);
77f951a8fcSStefano Zampini     corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs);
78f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,4,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
79454e267fSLisandro Dalcin   }
80454e267fSLisandro Dalcin   *nel = da->ne;
81454e267fSLisandro Dalcin   *nen = nn;
82454e267fSLisandro Dalcin   *e   = da->e;
83454e267fSLisandro Dalcin   PetscFunctionReturn(0);
84454e267fSLisandro Dalcin }
85454e267fSLisandro Dalcin 
862dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
87454e267fSLisandro Dalcin {
88454e267fSLisandro Dalcin   PetscErrorCode ierr;
89454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
90454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
91454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
92454e267fSLisandro Dalcin   PetscInt       k,zs,ze,Zs,Ze;
93454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[8], ns=6, nn=4;
94454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
95454e267fSLisandro Dalcin                                0,1,7,4,
96454e267fSLisandro Dalcin                                1,2,3,7,
97454e267fSLisandro Dalcin                                1,2,7,6,
98454e267fSLisandro Dalcin                                1,4,5,7,
99454e267fSLisandro Dalcin                                1,5,6,7};
1005fd66863SKarl Rupp 
101454e267fSLisandro Dalcin   PetscFunctionBegin;
102b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_P1) {nn=4;}
103b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_Q1) {nn=8;}
104454e267fSLisandro Dalcin   if (!da->e) {
105f951a8fcSStefano Zampini     PetscInt corners[8];
106f951a8fcSStefano Zampini 
10798b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
108b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6;}
109b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
110454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
111454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
112454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
113454e267fSLisandro Dalcin     ye    += ys; Ye += Ys; if (ys != Ys) ys -= 1;
114454e267fSLisandro Dalcin     ze    += zs; Ze += Zs; if (zs != Zs) zs -= 1;
115454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
116854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr);
117454e267fSLisandro Dalcin     for (k=zs; k<ze-1; k++) {
118454e267fSLisandro Dalcin       for (j=ys; j<ye-1; j++) {
119454e267fSLisandro Dalcin         for (i=xs; i<xe-1; i++) {
120454e267fSLisandro Dalcin           cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
121454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
122454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
123454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
124454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
125454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
126454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
127454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
128454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
1298865f1eaSKarl Rupp             for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]];
130454e267fSLisandro Dalcin           }
131454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
1328865f1eaSKarl Rupp             for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c];
133454e267fSLisandro Dalcin           }
134454e267fSLisandro Dalcin         }
135454e267fSLisandro Dalcin       }
136454e267fSLisandro Dalcin     }
137f951a8fcSStefano Zampini     corners[0] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
138f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
139f951a8fcSStefano Zampini     corners[2] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
140f951a8fcSStefano Zampini     corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
141f951a8fcSStefano Zampini     corners[4] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
142f951a8fcSStefano Zampini     corners[5] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
143f951a8fcSStefano Zampini     corners[6] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
144f951a8fcSStefano Zampini     corners[7] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
145f951a8fcSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,8,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr);
146454e267fSLisandro Dalcin   }
147454e267fSLisandro Dalcin   *nel = da->ne;
148454e267fSLisandro Dalcin   *nen = nn;
149454e267fSLisandro Dalcin   *e   = da->e;
150454e267fSLisandro Dalcin   PetscFunctionReturn(0);
151454e267fSLisandro Dalcin }
152454e267fSLisandro Dalcin 
153f5f57ec0SBarry Smith /*@
154d4a6ed37SStefano Zampini    DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left
155d4a6ed37SStefano Zampini    corner of the non-overlapping decomposition identified by DMDAGetElements()
156d4a6ed37SStefano Zampini 
157d4a6ed37SStefano Zampini     Not Collective
158d4a6ed37SStefano Zampini 
159d4a6ed37SStefano Zampini    Input Parameter:
160d4a6ed37SStefano Zampini .     da - the DM object
161d4a6ed37SStefano Zampini 
162d4a6ed37SStefano Zampini    Output Parameters:
163d4a6ed37SStefano Zampini +     gx - the x index
164d4a6ed37SStefano Zampini .     gy - the y index
165d4a6ed37SStefano Zampini -     gz - the z index
166d4a6ed37SStefano Zampini 
167d4a6ed37SStefano Zampini    Level: intermediate
168d4a6ed37SStefano Zampini 
169d4a6ed37SStefano Zampini    Notes:
170d4a6ed37SStefano Zampini 
171d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
172d4a6ed37SStefano Zampini @*/
173d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz)
174d4a6ed37SStefano Zampini {
175d4a6ed37SStefano Zampini   PetscInt       xs,Xs;
176d4a6ed37SStefano Zampini   PetscInt       ys,Ys;
177d4a6ed37SStefano Zampini   PetscInt       zs,Zs;
178d4a6ed37SStefano Zampini   PetscBool      isda;
179d4a6ed37SStefano Zampini   PetscErrorCode ierr;
180d4a6ed37SStefano Zampini 
181d4a6ed37SStefano Zampini   PetscFunctionBegin;
182d4a6ed37SStefano Zampini   PetscValidHeaderSpecific(da,DM_CLASSID,1);
183d4a6ed37SStefano Zampini   if (gx) PetscValidIntPointer(gx,2);
184d4a6ed37SStefano Zampini   if (gy) PetscValidIntPointer(gy,3);
185d4a6ed37SStefano Zampini   if (gz) PetscValidIntPointer(gz,4);
186d4a6ed37SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
187d4a6ed37SStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
188d4a6ed37SStefano Zampini   ierr = DMDAGetCorners(da,&xs,&ys,&zs,NULL,NULL,NULL);CHKERRQ(ierr);
189d4a6ed37SStefano Zampini   ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr);
190d4a6ed37SStefano Zampini   if (xs != Xs) xs -= 1;
191d4a6ed37SStefano Zampini   if (ys != Ys) ys -= 1;
192d4a6ed37SStefano Zampini   if (zs != Zs) zs -= 1;
193d4a6ed37SStefano Zampini   if (gx) *gx  = xs;
194d4a6ed37SStefano Zampini   if (gy) *gy  = ys;
195d4a6ed37SStefano Zampini   if (gz) *gz  = zs;
196d4a6ed37SStefano Zampini   PetscFunctionReturn(0);
197d4a6ed37SStefano Zampini }
198d4a6ed37SStefano Zampini 
199d4a6ed37SStefano Zampini /*@
200d4a6ed37SStefano Zampini       DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by DMDAGetElements()
201d4a6ed37SStefano Zampini 
202d4a6ed37SStefano Zampini     Not Collective
203d4a6ed37SStefano Zampini 
204d4a6ed37SStefano Zampini    Input Parameter:
205d4a6ed37SStefano Zampini .     da - the DM object
206d4a6ed37SStefano Zampini 
207d4a6ed37SStefano Zampini    Output Parameters:
208d4a6ed37SStefano Zampini +     mx - number of local elements in x-direction
209d4a6ed37SStefano Zampini .     my - number of local elements in y-direction
210d4a6ed37SStefano Zampini -     mz - number of local elements in z-direction
211d4a6ed37SStefano Zampini 
212d4a6ed37SStefano Zampini    Level: intermediate
213d4a6ed37SStefano Zampini 
214*95452b02SPatrick Sanan    Notes:
215*95452b02SPatrick Sanan     It returns the same number of elements, irrespective of the DMDAElementType
216d4a6ed37SStefano Zampini 
217d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements
218d4a6ed37SStefano Zampini @*/
219d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz)
220d4a6ed37SStefano Zampini {
221d4a6ed37SStefano Zampini   PetscInt       xs,xe,Xs;
222d4a6ed37SStefano Zampini   PetscInt       ys,ye,Ys;
223d4a6ed37SStefano Zampini   PetscInt       zs,ze,Zs;
224d4a6ed37SStefano Zampini   PetscInt       dim;
225d4a6ed37SStefano Zampini   PetscBool      isda;
226d4a6ed37SStefano Zampini   PetscErrorCode ierr;
227d4a6ed37SStefano Zampini 
228d4a6ed37SStefano Zampini   PetscFunctionBegin;
229d4a6ed37SStefano Zampini   PetscValidHeaderSpecific(da,DM_CLASSID,1);
230d4a6ed37SStefano Zampini   if (mx) PetscValidIntPointer(mx,2);
231d4a6ed37SStefano Zampini   if (my) PetscValidIntPointer(my,3);
232d4a6ed37SStefano Zampini   if (mz) PetscValidIntPointer(mz,4);
233d4a6ed37SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
234d4a6ed37SStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
235d4a6ed37SStefano Zampini   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
236d4a6ed37SStefano Zampini   ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr);
237d4a6ed37SStefano Zampini   xe  += xs; if (xs != Xs) xs -= 1;
238d4a6ed37SStefano Zampini   ye  += ys; if (ys != Ys) ys -= 1;
239d4a6ed37SStefano Zampini   ze  += zs; if (zs != Zs) zs -= 1;
240d4a6ed37SStefano Zampini   if (mx) *mx  = 0;
241d4a6ed37SStefano Zampini   if (my) *my  = 0;
242d4a6ed37SStefano Zampini   if (mz) *mz  = 0;
243d4a6ed37SStefano Zampini   ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
244d4a6ed37SStefano Zampini   switch (dim) {
245d4a6ed37SStefano Zampini   case 3:
246d4a6ed37SStefano Zampini     if (mz) *mz = ze - zs - 1;
247d4a6ed37SStefano Zampini   case 2:
248d4a6ed37SStefano Zampini     if (my) *my = ye - ys - 1;
249d4a6ed37SStefano Zampini   case 1:
250d4a6ed37SStefano Zampini     if (mx) *mx = xe - xs - 1;
251d4a6ed37SStefano Zampini     break;
252d4a6ed37SStefano Zampini   }
253d4a6ed37SStefano Zampini   PetscFunctionReturn(0);
254d4a6ed37SStefano Zampini }
255d4a6ed37SStefano Zampini 
256d4a6ed37SStefano Zampini /*@
2572dde6fd4SLisandro Dalcin       DMDASetElementType - Sets the element type to be returned by DMDAGetElements()
2582dde6fd4SLisandro Dalcin 
2592dde6fd4SLisandro Dalcin     Not Collective
2602dde6fd4SLisandro Dalcin 
2612dde6fd4SLisandro Dalcin    Input Parameter:
2622dde6fd4SLisandro Dalcin .     da - the DMDA object
2632dde6fd4SLisandro Dalcin 
2642dde6fd4SLisandro Dalcin    Output Parameters:
2656420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
2662dde6fd4SLisandro Dalcin 
2672dde6fd4SLisandro Dalcin    Level: intermediate
2682dde6fd4SLisandro Dalcin 
2692dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
2702dde6fd4SLisandro Dalcin @*/
2712dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
2722dde6fd4SLisandro Dalcin {
2732dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
2742dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
275f951a8fcSStefano Zampini   PetscBool      isda;
2762dde6fd4SLisandro Dalcin 
2772dde6fd4SLisandro Dalcin   PetscFunctionBegin;
2782dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2792dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
280f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
281f951a8fcSStefano Zampini   if (!isda) PetscFunctionReturn(0);
2822dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
2832dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
284f951a8fcSStefano Zampini     ierr = ISDestroy(&dd->ecorners);CHKERRQ(ierr);
2858865f1eaSKarl Rupp 
2862dde6fd4SLisandro Dalcin     dd->elementtype = etype;
2872dde6fd4SLisandro Dalcin     dd->ne          = 0;
2880298fd71SBarry Smith     dd->e           = NULL;
2892dde6fd4SLisandro Dalcin   }
2902dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
2912dde6fd4SLisandro Dalcin }
2922dde6fd4SLisandro Dalcin 
293f5f57ec0SBarry Smith /*@
2942dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
2952dde6fd4SLisandro Dalcin 
2962dde6fd4SLisandro Dalcin     Not Collective
2972dde6fd4SLisandro Dalcin 
2982dde6fd4SLisandro Dalcin    Input Parameter:
2992dde6fd4SLisandro Dalcin .     da - the DMDA object
3002dde6fd4SLisandro Dalcin 
3012dde6fd4SLisandro Dalcin    Output Parameters:
3026420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
3032dde6fd4SLisandro Dalcin 
3042dde6fd4SLisandro Dalcin    Level: intermediate
3052dde6fd4SLisandro Dalcin 
3062dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
3072dde6fd4SLisandro Dalcin @*/
3082dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
3092dde6fd4SLisandro Dalcin {
3102dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
311f951a8fcSStefano Zampini   PetscErrorCode ierr;
312f951a8fcSStefano Zampini   PetscBool      isda;
3135fd66863SKarl Rupp 
3142dde6fd4SLisandro Dalcin   PetscFunctionBegin;
3152dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
3162dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
317f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
318f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
3192dde6fd4SLisandro Dalcin   *etype = dd->elementtype;
3202dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
3212dde6fd4SLisandro Dalcin }
3222dde6fd4SLisandro Dalcin 
3232dde6fd4SLisandro Dalcin /*@C
3242dde6fd4SLisandro Dalcin       DMDAGetElements - Gets an array containing the indices (in local coordinates)
3252dde6fd4SLisandro Dalcin                  of all the local elements
3262dde6fd4SLisandro Dalcin 
3272dde6fd4SLisandro Dalcin     Not Collective
3282dde6fd4SLisandro Dalcin 
3292dde6fd4SLisandro Dalcin    Input Parameter:
3302dde6fd4SLisandro Dalcin .     dm - the DM object
3312dde6fd4SLisandro Dalcin 
3322dde6fd4SLisandro Dalcin    Output Parameters:
3332dde6fd4SLisandro Dalcin +     nel - number of local elements
3342dde6fd4SLisandro Dalcin .     nen - number of element nodes
335c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
3362dde6fd4SLisandro Dalcin 
3372dde6fd4SLisandro Dalcin    Level: intermediate
3382dde6fd4SLisandro Dalcin 
33993386343SJed Brown    Notes:
34093386343SJed Brown      Call DMDARestoreElements() once you have finished accessing the elements.
34193386343SJed Brown 
34245950774SSatish Balay      Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.
343009f0576SBarry Smith 
344009f0576SBarry 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.
345009f0576SBarry Smith 
346f5f57ec0SBarry Smith      Not supported in Fortran
347f5f57ec0SBarry Smith 
34893386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
3492dde6fd4SLisandro Dalcin @*/
3502dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
351454e267fSLisandro Dalcin {
352c73cfb54SMatthew G. Knepley   PetscInt       dim;
353454e267fSLisandro Dalcin   PetscErrorCode ierr;
35499c57625SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
355f951a8fcSStefano Zampini   PetscBool      isda;
3565fd66863SKarl Rupp 
357454e267fSLisandro Dalcin   PetscFunctionBegin;
358f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
359f951a8fcSStefano Zampini   PetscValidIntPointer(nel,2);
360f951a8fcSStefano Zampini   PetscValidIntPointer(nen,3);
361f951a8fcSStefano Zampini   PetscValidPointer(e,4);
362f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
363f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
364f951a8fcSStefano 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");
365c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
366c73cfb54SMatthew G. Knepley   if (dim==-1) {
3670298fd71SBarry Smith     *nel = 0; *nen = 0; *e = NULL;
368c73cfb54SMatthew G. Knepley   } else if (dim==1) {
3692dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
370c73cfb54SMatthew G. Knepley   } else if (dim==2) {
3712dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
372c73cfb54SMatthew G. Knepley   } else if (dim==3) {
3732dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
374c73cfb54SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
375454e267fSLisandro Dalcin   PetscFunctionReturn(0);
376454e267fSLisandro Dalcin }
3772dde6fd4SLisandro Dalcin 
378f951a8fcSStefano Zampini /*@
379d4a6ed37SStefano Zampini       DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates)
380f951a8fcSStefano Zampini                                  of the non-overlapping decomposition identified by DMDAGetElements
381f951a8fcSStefano Zampini 
382f951a8fcSStefano Zampini     Not Collective
383f951a8fcSStefano Zampini 
384f951a8fcSStefano Zampini    Input Parameter:
385f951a8fcSStefano Zampini .     dm - the DM object
386f951a8fcSStefano Zampini 
387f951a8fcSStefano Zampini    Output Parameters:
388f951a8fcSStefano Zampini .     is - the index set
389f951a8fcSStefano Zampini 
390f951a8fcSStefano Zampini    Level: intermediate
391f951a8fcSStefano Zampini 
392*95452b02SPatrick Sanan    Notes:
393*95452b02SPatrick Sanan     Call DMDARestoreSubdomainCornersIS() once you have finished accessing the index set.
394f951a8fcSStefano Zampini 
395f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElementsCornersIS()
396f951a8fcSStefano Zampini @*/
397d4a6ed37SStefano Zampini PetscErrorCode  DMDAGetSubdomainCornersIS(DM dm,IS *is)
398f951a8fcSStefano Zampini {
399f951a8fcSStefano Zampini   PetscErrorCode ierr;
400f951a8fcSStefano Zampini   DM_DA          *dd = (DM_DA*)dm->data;
401f951a8fcSStefano Zampini   PetscBool      isda;
402f951a8fcSStefano Zampini 
403f951a8fcSStefano Zampini   PetscFunctionBegin;
404f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
405f951a8fcSStefano Zampini   PetscValidPointer(is,2);
406f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
407f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
408f951a8fcSStefano 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");
409f951a8fcSStefano Zampini   if (!dd->ecorners) { /* compute elements if not yet done */
410f951a8fcSStefano Zampini     const PetscInt *e;
411f951a8fcSStefano Zampini     PetscInt       nel,nen;
412f951a8fcSStefano Zampini 
413f951a8fcSStefano Zampini     ierr = DMDAGetElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
414f951a8fcSStefano Zampini     ierr = DMDARestoreElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
415f951a8fcSStefano Zampini   }
416f951a8fcSStefano Zampini   *is = dd->ecorners;
417f951a8fcSStefano Zampini   PetscFunctionReturn(0);
418f951a8fcSStefano Zampini }
419f951a8fcSStefano Zampini 
4202dde6fd4SLisandro Dalcin /*@C
421009f0576SBarry Smith       DMDARestoreElements - Restores the array obtained with DMDAGetElements()
4222dde6fd4SLisandro Dalcin 
4232dde6fd4SLisandro Dalcin     Not Collective
4242dde6fd4SLisandro Dalcin 
4252dde6fd4SLisandro Dalcin    Input Parameter:
4262dde6fd4SLisandro Dalcin +     dm - the DM object
4272dde6fd4SLisandro Dalcin .     nel - number of local elements
4282dde6fd4SLisandro Dalcin .     nen - number of element nodes
429c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
4302dde6fd4SLisandro Dalcin 
4312dde6fd4SLisandro Dalcin    Level: intermediate
4322dde6fd4SLisandro Dalcin 
433009f0576SBarry Smith    Note: You should not access these values after you have called this routine.
434009f0576SBarry Smith 
435c67aec63SBarry Smith          This restore signals the DMDA object that you no longer need access to the array information.
436009f0576SBarry Smith 
437f5f57ec0SBarry Smith          Not supported in Fortran
438f5f57ec0SBarry Smith 
4392dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
4402dde6fd4SLisandro Dalcin @*/
4412dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
4422dde6fd4SLisandro Dalcin {
4432dde6fd4SLisandro Dalcin   PetscFunctionBegin;
4442dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4452dde6fd4SLisandro Dalcin   PetscValidIntPointer(nel,2);
4462dde6fd4SLisandro Dalcin   PetscValidIntPointer(nen,3);
4472dde6fd4SLisandro Dalcin   PetscValidPointer(e,4);
4489cc8bc54SJed Brown   *nel = 0;
4499cc8bc54SJed Brown   *nen = -1;
4509cc8bc54SJed Brown   *e = NULL;
4512dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
4522dde6fd4SLisandro Dalcin }
453f951a8fcSStefano Zampini 
454f951a8fcSStefano Zampini /*@
455d4a6ed37SStefano Zampini       DMDARestoreSubdomainCornersIS - Restores the IS obtained with DMDAGetSubdomainCornersIS()
456f951a8fcSStefano Zampini 
457f951a8fcSStefano Zampini     Not Collective
458f951a8fcSStefano Zampini 
459f951a8fcSStefano Zampini    Input Parameter:
460f951a8fcSStefano Zampini +     dm - the DM object
461f951a8fcSStefano Zampini -     is - the index set
462f951a8fcSStefano Zampini 
463f951a8fcSStefano Zampini    Level: intermediate
464f951a8fcSStefano Zampini 
465f951a8fcSStefano Zampini    Note:
466f951a8fcSStefano Zampini 
467d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetSubdomainCornersIS()
468f951a8fcSStefano Zampini @*/
469d4a6ed37SStefano Zampini PetscErrorCode  DMDARestoreSubdomainCornersIS(DM dm,IS *is)
470f951a8fcSStefano Zampini {
471f951a8fcSStefano Zampini   PetscFunctionBegin;
472f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
473f951a8fcSStefano Zampini   PetscValidHeaderSpecific(*is,IS_CLASSID,2);
474f951a8fcSStefano Zampini   *is = NULL;
475f951a8fcSStefano Zampini   PetscFunctionReturn(0);
476f951a8fcSStefano Zampini }
477