xref: /petsc/src/dm/impls/da/dagetelem.c (revision f951a8fcfcf7163ef637f761a5651bb04e634f4d)
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) {
13*f951a8fcSStefano Zampini     PetscInt corners[2];
14*f951a8fcSStefano 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*f951a8fcSStefano Zampini     corners[0] = (xs  -Xs);
26*f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs);
27*f951a8fcSStefano 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) {
49*f951a8fcSStefano Zampini     PetscInt corners[4];
50*f951a8fcSStefano 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     }
74*f951a8fcSStefano Zampini     corners[0] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs);
75*f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs);
76*f951a8fcSStefano Zampini     corners[2] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs);
77*f951a8fcSStefano Zampini     corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs);
78*f951a8fcSStefano 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) {
105*f951a8fcSStefano Zampini     PetscInt corners[8];
106*f951a8fcSStefano 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     }
137*f951a8fcSStefano Zampini     corners[0] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
138*f951a8fcSStefano Zampini     corners[1] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
139*f951a8fcSStefano Zampini     corners[2] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
140*f951a8fcSStefano Zampini     corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (zs-Zs  )*(Xe-Xs)*(Ye-Ys);
141*f951a8fcSStefano Zampini     corners[4] = (xs  -Xs) + (ys  -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
142*f951a8fcSStefano Zampini     corners[5] = (xe-1-Xs) + (ys  -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
143*f951a8fcSStefano Zampini     corners[6] = (xs  -Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
144*f951a8fcSStefano Zampini     corners[7] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys);
145*f951a8fcSStefano 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 /*@
1542dde6fd4SLisandro Dalcin       DMDASetElementType - Sets the element type to be returned by DMDAGetElements()
1552dde6fd4SLisandro Dalcin 
1562dde6fd4SLisandro Dalcin     Not Collective
1572dde6fd4SLisandro Dalcin 
1582dde6fd4SLisandro Dalcin    Input Parameter:
1592dde6fd4SLisandro Dalcin .     da - the DMDA object
1602dde6fd4SLisandro Dalcin 
1612dde6fd4SLisandro Dalcin    Output Parameters:
1626420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
1632dde6fd4SLisandro Dalcin 
1642dde6fd4SLisandro Dalcin    Level: intermediate
1652dde6fd4SLisandro Dalcin 
1662dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
1672dde6fd4SLisandro Dalcin @*/
1682dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
1692dde6fd4SLisandro Dalcin {
1702dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
1712dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
172*f951a8fcSStefano Zampini   PetscBool      isda;
1732dde6fd4SLisandro Dalcin 
1742dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1752dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1762dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
177*f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
178*f951a8fcSStefano Zampini   if (!isda) PetscFunctionReturn(0);
1792dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
1802dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
181*f951a8fcSStefano Zampini     ierr = ISDestroy(&dd->ecorners);CHKERRQ(ierr);
1828865f1eaSKarl Rupp 
1832dde6fd4SLisandro Dalcin     dd->elementtype = etype;
1842dde6fd4SLisandro Dalcin     dd->ne          = 0;
1850298fd71SBarry Smith     dd->e           = NULL;
1862dde6fd4SLisandro Dalcin   }
1872dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1882dde6fd4SLisandro Dalcin }
1892dde6fd4SLisandro Dalcin 
190f5f57ec0SBarry Smith /*@
1912dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
1922dde6fd4SLisandro Dalcin 
1932dde6fd4SLisandro Dalcin     Not Collective
1942dde6fd4SLisandro Dalcin 
1952dde6fd4SLisandro Dalcin    Input Parameter:
1962dde6fd4SLisandro Dalcin .     da - the DMDA object
1972dde6fd4SLisandro Dalcin 
1982dde6fd4SLisandro Dalcin    Output Parameters:
1996420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
2002dde6fd4SLisandro Dalcin 
2012dde6fd4SLisandro Dalcin    Level: intermediate
2022dde6fd4SLisandro Dalcin 
2032dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
2042dde6fd4SLisandro Dalcin @*/
2052dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
2062dde6fd4SLisandro Dalcin {
2072dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
208*f951a8fcSStefano Zampini   PetscErrorCode ierr;
209*f951a8fcSStefano Zampini   PetscBool      isda;
2105fd66863SKarl Rupp 
2112dde6fd4SLisandro Dalcin   PetscFunctionBegin;
2122dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2132dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
214*f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr);
215*f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name);
2162dde6fd4SLisandro Dalcin   *etype = dd->elementtype;
2172dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
2182dde6fd4SLisandro Dalcin }
2192dde6fd4SLisandro Dalcin 
2202dde6fd4SLisandro Dalcin /*@C
2212dde6fd4SLisandro Dalcin       DMDAGetElements - Gets an array containing the indices (in local coordinates)
2222dde6fd4SLisandro Dalcin                  of all the local elements
2232dde6fd4SLisandro Dalcin 
2242dde6fd4SLisandro Dalcin     Not Collective
2252dde6fd4SLisandro Dalcin 
2262dde6fd4SLisandro Dalcin    Input Parameter:
2272dde6fd4SLisandro Dalcin .     dm - the DM object
2282dde6fd4SLisandro Dalcin 
2292dde6fd4SLisandro Dalcin    Output Parameters:
2302dde6fd4SLisandro Dalcin +     nel - number of local elements
2312dde6fd4SLisandro Dalcin .     nen - number of element nodes
232c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2332dde6fd4SLisandro Dalcin 
2342dde6fd4SLisandro Dalcin    Level: intermediate
2352dde6fd4SLisandro Dalcin 
23693386343SJed Brown    Notes:
23793386343SJed Brown      Call DMDARestoreElements() once you have finished accessing the elements.
23893386343SJed Brown 
23945950774SSatish Balay      Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.
240009f0576SBarry Smith 
241009f0576SBarry 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.
242009f0576SBarry Smith 
243f5f57ec0SBarry Smith      Not supported in Fortran
244f5f57ec0SBarry Smith 
24593386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
2462dde6fd4SLisandro Dalcin @*/
2472dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
248454e267fSLisandro Dalcin {
249c73cfb54SMatthew G. Knepley   PetscInt       dim;
250454e267fSLisandro Dalcin   PetscErrorCode ierr;
25199c57625SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
252*f951a8fcSStefano Zampini   PetscBool      isda;
2535fd66863SKarl Rupp 
254454e267fSLisandro Dalcin   PetscFunctionBegin;
255*f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
256*f951a8fcSStefano Zampini   PetscValidIntPointer(nel,2);
257*f951a8fcSStefano Zampini   PetscValidIntPointer(nen,3);
258*f951a8fcSStefano Zampini   PetscValidPointer(e,4);
259*f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
260*f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
261*f951a8fcSStefano 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");
262c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
263c73cfb54SMatthew G. Knepley   if (dim==-1) {
2640298fd71SBarry Smith     *nel = 0; *nen = 0; *e = NULL;
265c73cfb54SMatthew G. Knepley   } else if (dim==1) {
2662dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
267c73cfb54SMatthew G. Knepley   } else if (dim==2) {
2682dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
269c73cfb54SMatthew G. Knepley   } else if (dim==3) {
2702dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
271c73cfb54SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
272454e267fSLisandro Dalcin   PetscFunctionReturn(0);
273454e267fSLisandro Dalcin }
2742dde6fd4SLisandro Dalcin 
275*f951a8fcSStefano Zampini /*@
276*f951a8fcSStefano Zampini       DMDAGetElementsCornersIS - Gets an index set containing the corner indices (in local coordinates)
277*f951a8fcSStefano Zampini                                  of the non-overlapping decomposition identified by DMDAGetElements
278*f951a8fcSStefano Zampini 
279*f951a8fcSStefano Zampini     Not Collective
280*f951a8fcSStefano Zampini 
281*f951a8fcSStefano Zampini    Input Parameter:
282*f951a8fcSStefano Zampini .     dm - the DM object
283*f951a8fcSStefano Zampini 
284*f951a8fcSStefano Zampini    Output Parameters:
285*f951a8fcSStefano Zampini .     is - the index set
286*f951a8fcSStefano Zampini 
287*f951a8fcSStefano Zampini    Level: intermediate
288*f951a8fcSStefano Zampini 
289*f951a8fcSStefano Zampini    Notes: Call DMDARestoreElementsCornersIS() once you have finished accessing the index set.
290*f951a8fcSStefano Zampini 
291*f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElementsCornersIS()
292*f951a8fcSStefano Zampini @*/
293*f951a8fcSStefano Zampini PetscErrorCode  DMDAGetElementsCornersIS(DM dm,IS *is)
294*f951a8fcSStefano Zampini {
295*f951a8fcSStefano Zampini   PetscErrorCode ierr;
296*f951a8fcSStefano Zampini   DM_DA          *dd = (DM_DA*)dm->data;
297*f951a8fcSStefano Zampini   PetscBool      isda;
298*f951a8fcSStefano Zampini 
299*f951a8fcSStefano Zampini   PetscFunctionBegin;
300*f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
301*f951a8fcSStefano Zampini   PetscValidPointer(is,2);
302*f951a8fcSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr);
303*f951a8fcSStefano Zampini   if (!isda) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name);
304*f951a8fcSStefano 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");
305*f951a8fcSStefano Zampini   if (!dd->ecorners) { /* compute elements if not yet done */
306*f951a8fcSStefano Zampini     const PetscInt *e;
307*f951a8fcSStefano Zampini     PetscInt       nel,nen;
308*f951a8fcSStefano Zampini 
309*f951a8fcSStefano Zampini     ierr = DMDAGetElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
310*f951a8fcSStefano Zampini     ierr = DMDARestoreElements(dm,&nel,&nen,&e);CHKERRQ(ierr);
311*f951a8fcSStefano Zampini   }
312*f951a8fcSStefano Zampini   *is = dd->ecorners;
313*f951a8fcSStefano Zampini   PetscFunctionReturn(0);
314*f951a8fcSStefano Zampini }
315*f951a8fcSStefano Zampini 
3162dde6fd4SLisandro Dalcin /*@C
317009f0576SBarry Smith       DMDARestoreElements - Restores the array obtained with DMDAGetElements()
3182dde6fd4SLisandro Dalcin 
3192dde6fd4SLisandro Dalcin     Not Collective
3202dde6fd4SLisandro Dalcin 
3212dde6fd4SLisandro Dalcin    Input Parameter:
3222dde6fd4SLisandro Dalcin +     dm - the DM object
3232dde6fd4SLisandro Dalcin .     nel - number of local elements
3242dde6fd4SLisandro Dalcin .     nen - number of element nodes
325c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
3262dde6fd4SLisandro Dalcin 
3272dde6fd4SLisandro Dalcin    Level: intermediate
3282dde6fd4SLisandro Dalcin 
329009f0576SBarry Smith    Note: You should not access these values after you have called this routine.
330009f0576SBarry Smith 
331c67aec63SBarry Smith          This restore signals the DMDA object that you no longer need access to the array information.
332009f0576SBarry Smith 
333f5f57ec0SBarry Smith          Not supported in Fortran
334f5f57ec0SBarry Smith 
3352dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
3362dde6fd4SLisandro Dalcin @*/
3372dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
3382dde6fd4SLisandro Dalcin {
3392dde6fd4SLisandro Dalcin   PetscFunctionBegin;
3402dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3412dde6fd4SLisandro Dalcin   PetscValidIntPointer(nel,2);
3422dde6fd4SLisandro Dalcin   PetscValidIntPointer(nen,3);
3432dde6fd4SLisandro Dalcin   PetscValidPointer(e,4);
3449cc8bc54SJed Brown   *nel = 0;
3459cc8bc54SJed Brown   *nen = -1;
3469cc8bc54SJed Brown   *e = NULL;
3472dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
3482dde6fd4SLisandro Dalcin }
349*f951a8fcSStefano Zampini 
350*f951a8fcSStefano Zampini /*@
351*f951a8fcSStefano Zampini       DMDARestoreElementsCornersIS - Restores the IS obtained with DMDAGetElementsCornersIS()
352*f951a8fcSStefano Zampini 
353*f951a8fcSStefano Zampini     Not Collective
354*f951a8fcSStefano Zampini 
355*f951a8fcSStefano Zampini    Input Parameter:
356*f951a8fcSStefano Zampini +     dm - the DM object
357*f951a8fcSStefano Zampini -     is - the index set
358*f951a8fcSStefano Zampini 
359*f951a8fcSStefano Zampini    Level: intermediate
360*f951a8fcSStefano Zampini 
361*f951a8fcSStefano Zampini    Note:
362*f951a8fcSStefano Zampini 
363*f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElementsCornersIS()
364*f951a8fcSStefano Zampini @*/
365*f951a8fcSStefano Zampini PetscErrorCode  DMDARestoreElementsCornersIS(DM dm,IS *is)
366*f951a8fcSStefano Zampini {
367*f951a8fcSStefano Zampini   PetscFunctionBegin;
368*f951a8fcSStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
369*f951a8fcSStefano Zampini   PetscValidHeaderSpecific(*is,IS_CLASSID,2);
370*f951a8fcSStefano Zampini   *is = NULL;
371*f951a8fcSStefano Zampini   PetscFunctionReturn(0);
372*f951a8fcSStefano Zampini }
373