xref: /petsc/src/dm/impls/da/dagetelem.c (revision c2cd2aa3969a14fa686ea977e3dce69115efdab6)
1454e267fSLisandro Dalcin 
2c6db04a5SJed Brown #include <private/daimpl.h>     /*I  "petscdmda.h"   I*/
3454e267fSLisandro Dalcin 
4454e267fSLisandro Dalcin #undef __FUNCT__
52dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_1D"
62dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_1D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
7454e267fSLisandro Dalcin {
8454e267fSLisandro Dalcin   PetscErrorCode ierr;
9454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
10454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
11454e267fSLisandro Dalcin   PetscInt       cnt=0;
12454e267fSLisandro Dalcin   PetscFunctionBegin;
13454e267fSLisandro Dalcin   if (!da->e) {
14454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,0,0,&xe,0,0);CHKERRQ(ierr);
15454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,0,0,&Xe,0,0);CHKERRQ(ierr);
16454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
17454e267fSLisandro Dalcin     da->ne = 1*(xe - xs - 1);
18454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + 2*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
19454e267fSLisandro Dalcin     for (i=xs; i<xe-1; i++) {
20454e267fSLisandro Dalcin       da->e[cnt++] = (i-Xs  );
21454e267fSLisandro Dalcin       da->e[cnt++] = (i-Xs+1);
22454e267fSLisandro Dalcin     }
23454e267fSLisandro Dalcin   }
24454e267fSLisandro Dalcin   *nel = da->ne;
25454e267fSLisandro Dalcin   *nen = 2;
26454e267fSLisandro Dalcin   *e   = da->e;
27454e267fSLisandro Dalcin   PetscFunctionReturn(0);
28454e267fSLisandro Dalcin }
29454e267fSLisandro Dalcin 
30454e267fSLisandro Dalcin #undef __FUNCT__
312dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_2D"
322dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
33454e267fSLisandro Dalcin {
34454e267fSLisandro Dalcin   PetscErrorCode ierr;
35454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
36454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
37454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
38454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[4], ns=2, nn=3;
39454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,
40454e267fSLisandro Dalcin                                2,3,1};
41454e267fSLisandro Dalcin   PetscFunctionBegin;
42454e267fSLisandro Dalcin   if (!da->e) {
43454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2; nn=3;}
44454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=4;}
45454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr);
46454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr);
47454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
48454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
49454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1);
50454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
51454e267fSLisandro Dalcin     for (j=ys; j<ye-1; j++) {
52454e267fSLisandro Dalcin       for (i=xs; i<xe-1; i++) {
53454e267fSLisandro Dalcin         cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs);
54454e267fSLisandro Dalcin         cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs);
55454e267fSLisandro Dalcin         cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs);
56454e267fSLisandro Dalcin         cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs);
57454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
58454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
59454e267fSLisandro Dalcin             da->e[cnt++] = cell[split[c]];
60454e267fSLisandro Dalcin         }
61454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
62454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
63454e267fSLisandro Dalcin             da->e[cnt++] = cell[c];
64454e267fSLisandro Dalcin         }
65454e267fSLisandro Dalcin       }
66454e267fSLisandro Dalcin     }
67454e267fSLisandro Dalcin   }
68454e267fSLisandro Dalcin   *nel = da->ne;
69454e267fSLisandro Dalcin   *nen = nn;
70454e267fSLisandro Dalcin   *e   = da->e;
71454e267fSLisandro Dalcin   PetscFunctionReturn(0);
72454e267fSLisandro Dalcin }
73454e267fSLisandro Dalcin 
74454e267fSLisandro Dalcin #undef __FUNCT__
752dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_3D"
762dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
77454e267fSLisandro Dalcin {
78454e267fSLisandro Dalcin   PetscErrorCode ierr;
79454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
80454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
81454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
82454e267fSLisandro Dalcin   PetscInt       k,zs,ze,Zs,Ze;
83454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[8], ns=6, nn=4;
84454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
85454e267fSLisandro Dalcin                                0,1,7,4,
86454e267fSLisandro Dalcin                                1,2,3,7,
87454e267fSLisandro Dalcin                                1,2,7,6,
88454e267fSLisandro Dalcin                                1,4,5,7,
89454e267fSLisandro Dalcin                                1,5,6,7};
90454e267fSLisandro Dalcin   PetscFunctionBegin;
91454e267fSLisandro Dalcin   if (!da->e) {
92454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6; nn=4;}
93454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=8;}
94454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
95454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
96454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
97454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
98454e267fSLisandro Dalcin     ze += zs; Ze += Zs; if (zs != Zs) zs -= 1;
99454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
100454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
101454e267fSLisandro Dalcin     for (k=zs; k<ze-1; k++) {
102454e267fSLisandro Dalcin       for (j=ys; j<ye-1; j++) {
103454e267fSLisandro Dalcin         for (i=xs; i<xe-1; i++) {
104454e267fSLisandro Dalcin           cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
105454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
106454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
107454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
108454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
109454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
110454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
111454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
112454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
113454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
114454e267fSLisandro Dalcin               da->e[cnt++] = cell[split[c]];
115454e267fSLisandro Dalcin           }
116454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
117454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
118454e267fSLisandro Dalcin               da->e[cnt++] = cell[c];
119454e267fSLisandro Dalcin           }
120454e267fSLisandro Dalcin         }
121454e267fSLisandro Dalcin       }
122454e267fSLisandro Dalcin     }
123454e267fSLisandro Dalcin   }
124454e267fSLisandro Dalcin   *nel = da->ne;
125454e267fSLisandro Dalcin   *nen = nn;
126454e267fSLisandro Dalcin   *e   = da->e;
127454e267fSLisandro Dalcin   PetscFunctionReturn(0);
128454e267fSLisandro Dalcin }
129454e267fSLisandro Dalcin 
1302dde6fd4SLisandro Dalcin /*@C
1312dde6fd4SLisandro Dalcin       DMDASetElementType - Sets the element type to be returned by DMDAGetElements()
1322dde6fd4SLisandro Dalcin 
1332dde6fd4SLisandro Dalcin     Not Collective
1342dde6fd4SLisandro Dalcin 
1352dde6fd4SLisandro Dalcin    Input Parameter:
1362dde6fd4SLisandro Dalcin .     da - the DMDA object
1372dde6fd4SLisandro Dalcin 
1382dde6fd4SLisandro Dalcin    Output Parameters:
1392dde6fd4SLisandro Dalcin .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
1402dde6fd4SLisandro Dalcin 
1412dde6fd4SLisandro Dalcin    Level: intermediate
1422dde6fd4SLisandro Dalcin 
1432dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
1442dde6fd4SLisandro Dalcin @*/
145454e267fSLisandro Dalcin #undef __FUNCT__
1462dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDASetElementType"
1472dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
1482dde6fd4SLisandro Dalcin {
1492dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
1502dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
1512dde6fd4SLisandro Dalcin 
1522dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1532dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1542dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
1552dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
1562dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
1572dde6fd4SLisandro Dalcin     dd->elementtype = etype;
1582dde6fd4SLisandro Dalcin     dd->ne          = 0;
1592dde6fd4SLisandro Dalcin     dd->e           = PETSC_NULL;
1602dde6fd4SLisandro Dalcin   }
1612dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1622dde6fd4SLisandro Dalcin }
1632dde6fd4SLisandro Dalcin 
1642dde6fd4SLisandro Dalcin /*@C
1652dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
1662dde6fd4SLisandro Dalcin 
1672dde6fd4SLisandro Dalcin     Not Collective
1682dde6fd4SLisandro Dalcin 
1692dde6fd4SLisandro Dalcin    Input Parameter:
1702dde6fd4SLisandro Dalcin .     da - the DMDA object
1712dde6fd4SLisandro Dalcin 
1722dde6fd4SLisandro Dalcin    Output Parameters:
1732dde6fd4SLisandro Dalcin .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
1742dde6fd4SLisandro Dalcin 
1752dde6fd4SLisandro Dalcin    Level: intermediate
1762dde6fd4SLisandro Dalcin 
1772dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
1782dde6fd4SLisandro Dalcin @*/
1792dde6fd4SLisandro Dalcin #undef __FUNCT__
1802dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElementType"
1812dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
1822dde6fd4SLisandro Dalcin {
1832dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
1842dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1852dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1862dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
1872dde6fd4SLisandro Dalcin   *etype = dd->elementtype;
1882dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1892dde6fd4SLisandro Dalcin }
1902dde6fd4SLisandro Dalcin 
1912dde6fd4SLisandro Dalcin /*@C
1922dde6fd4SLisandro Dalcin       DMDAGetElements - Gets an array containing the indices (in local coordinates)
1932dde6fd4SLisandro Dalcin                  of all the local elements
1942dde6fd4SLisandro Dalcin 
1952dde6fd4SLisandro Dalcin     Not Collective
1962dde6fd4SLisandro Dalcin 
1972dde6fd4SLisandro Dalcin    Input Parameter:
1982dde6fd4SLisandro Dalcin .     dm - the DM object
1992dde6fd4SLisandro Dalcin 
2002dde6fd4SLisandro Dalcin    Output Parameters:
2012dde6fd4SLisandro Dalcin +     nel - number of local elements
2022dde6fd4SLisandro Dalcin .     nen - number of element nodes
203*c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2042dde6fd4SLisandro Dalcin 
2052dde6fd4SLisandro Dalcin    Level: intermediate
2062dde6fd4SLisandro Dalcin 
207*c2cd2aa3SJed Brown .seealso: DMDAElementType, DMDASetElementType(), DMDARestoreElements(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
2082dde6fd4SLisandro Dalcin @*/
2092dde6fd4SLisandro Dalcin #undef __FUNCT__
2102dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements"
2112dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
212454e267fSLisandro Dalcin {
213454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
214454e267fSLisandro Dalcin   PetscErrorCode ierr;
215454e267fSLisandro Dalcin   PetscFunctionBegin;
216454e267fSLisandro Dalcin   if (da->dim==-1) {
217454e267fSLisandro Dalcin     *nel = 0; *nen = 0; *e = PETSC_NULL;
218454e267fSLisandro Dalcin   } else if (da->dim==1) {
2192dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
220454e267fSLisandro Dalcin   } else if (da->dim==2) {
2212dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
222454e267fSLisandro Dalcin   } else if (da->dim==3) {
2232dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
224454e267fSLisandro Dalcin   } else {
225454e267fSLisandro Dalcin     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",da->dim);
226454e267fSLisandro Dalcin   }
227454e267fSLisandro Dalcin 
228454e267fSLisandro Dalcin   PetscFunctionReturn(0);
229454e267fSLisandro Dalcin }
2302dde6fd4SLisandro Dalcin 
2312dde6fd4SLisandro Dalcin #undef __FUNCT__
2322dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDARestoreElements"
2332dde6fd4SLisandro Dalcin /*@C
2342dde6fd4SLisandro Dalcin       DMDARestoreElements - Returns an array containing the indices (in local coordinates)
2352dde6fd4SLisandro Dalcin                  of all the local elements obtained with DMDAGetElements()
2362dde6fd4SLisandro Dalcin 
2372dde6fd4SLisandro Dalcin     Not Collective
2382dde6fd4SLisandro Dalcin 
2392dde6fd4SLisandro Dalcin    Input Parameter:
2402dde6fd4SLisandro Dalcin +     dm - the DM object
2412dde6fd4SLisandro Dalcin .     nel - number of local elements
2422dde6fd4SLisandro Dalcin .     nen - number of element nodes
243*c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2442dde6fd4SLisandro Dalcin 
2452dde6fd4SLisandro Dalcin    Level: intermediate
2462dde6fd4SLisandro Dalcin 
2472dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
2482dde6fd4SLisandro Dalcin @*/
2492dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
2502dde6fd4SLisandro Dalcin {
2512dde6fd4SLisandro Dalcin   PetscFunctionBegin;
2522dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2532dde6fd4SLisandro Dalcin   PetscValidIntPointer(nel,2);
2542dde6fd4SLisandro Dalcin   PetscValidIntPointer(nen,3);
2552dde6fd4SLisandro Dalcin   PetscValidPointer(e,4);
2562dde6fd4SLisandro Dalcin   /* XXX */
2572dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
2582dde6fd4SLisandro Dalcin }
259