xref: /petsc/src/dm/impls/da/dagetelem.c (revision 5fd668637986a8d8518383a9159eebc368e1d5b4)
1454e267fSLisandro Dalcin 
2b45d2f2cSJed Brown #include <petsc-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;
12*5fd66863SKarl Rupp 
13454e267fSLisandro Dalcin   PetscFunctionBegin;
14454e267fSLisandro Dalcin   if (!da->e) {
15454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,0,0,&xe,0,0);CHKERRQ(ierr);
16454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,0,0,&Xe,0,0);CHKERRQ(ierr);
17454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
18454e267fSLisandro Dalcin     da->ne = 1*(xe - xs - 1);
19454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + 2*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
20454e267fSLisandro Dalcin     for (i=xs; i<xe-1; i++) {
21454e267fSLisandro Dalcin       da->e[cnt++] = (i-Xs  );
22454e267fSLisandro Dalcin       da->e[cnt++] = (i-Xs+1);
23454e267fSLisandro Dalcin     }
24454e267fSLisandro Dalcin   }
25454e267fSLisandro Dalcin   *nel = da->ne;
26454e267fSLisandro Dalcin   *nen = 2;
27454e267fSLisandro Dalcin   *e   = da->e;
28454e267fSLisandro Dalcin   PetscFunctionReturn(0);
29454e267fSLisandro Dalcin }
30454e267fSLisandro Dalcin 
31454e267fSLisandro Dalcin #undef __FUNCT__
322dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_2D"
332dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
34454e267fSLisandro Dalcin {
35454e267fSLisandro Dalcin   PetscErrorCode ierr;
36454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
37454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
38454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
39454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[4], ns=2, nn=3;
40454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,
41454e267fSLisandro Dalcin                                2,3,1};
42*5fd66863SKarl Rupp 
43454e267fSLisandro Dalcin   PetscFunctionBegin;
44454e267fSLisandro Dalcin   if (!da->e) {
45454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2; nn=3;}
46454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=4;}
47454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr);
48454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr);
49454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
50454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
51454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1);
52454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
53454e267fSLisandro Dalcin     for (j=ys; j<ye-1; j++) {
54454e267fSLisandro Dalcin       for (i=xs; i<xe-1; i++) {
55454e267fSLisandro Dalcin         cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs);
56454e267fSLisandro Dalcin         cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs);
57454e267fSLisandro Dalcin         cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs);
58454e267fSLisandro Dalcin         cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs);
59454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
60454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
61454e267fSLisandro Dalcin             da->e[cnt++] = cell[split[c]];
62454e267fSLisandro Dalcin         }
63454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
64454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
65454e267fSLisandro Dalcin             da->e[cnt++] = cell[c];
66454e267fSLisandro Dalcin         }
67454e267fSLisandro Dalcin       }
68454e267fSLisandro Dalcin     }
69454e267fSLisandro Dalcin   }
70454e267fSLisandro Dalcin   *nel = da->ne;
71454e267fSLisandro Dalcin   *nen = nn;
72454e267fSLisandro Dalcin   *e   = da->e;
73454e267fSLisandro Dalcin   PetscFunctionReturn(0);
74454e267fSLisandro Dalcin }
75454e267fSLisandro Dalcin 
76454e267fSLisandro Dalcin #undef __FUNCT__
772dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_3D"
782dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
79454e267fSLisandro Dalcin {
80454e267fSLisandro Dalcin   PetscErrorCode ierr;
81454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
82454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
83454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
84454e267fSLisandro Dalcin   PetscInt       k,zs,ze,Zs,Ze;
85454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[8], ns=6, nn=4;
86454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
87454e267fSLisandro Dalcin                                0,1,7,4,
88454e267fSLisandro Dalcin                                1,2,3,7,
89454e267fSLisandro Dalcin                                1,2,7,6,
90454e267fSLisandro Dalcin                                1,4,5,7,
91454e267fSLisandro Dalcin                                1,5,6,7};
92*5fd66863SKarl Rupp 
93454e267fSLisandro Dalcin   PetscFunctionBegin;
94454e267fSLisandro Dalcin   if (!da->e) {
95454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6; nn=4;}
96454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=8;}
97454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
98454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
99454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
100454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
101454e267fSLisandro Dalcin     ze += zs; Ze += Zs; if (zs != Zs) zs -= 1;
102454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
103454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
104454e267fSLisandro Dalcin     for (k=zs; k<ze-1; k++) {
105454e267fSLisandro Dalcin       for (j=ys; j<ye-1; j++) {
106454e267fSLisandro Dalcin         for (i=xs; i<xe-1; i++) {
107454e267fSLisandro Dalcin           cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
108454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
109454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
110454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
111454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
112454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
113454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
114454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
115454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
116454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
117454e267fSLisandro Dalcin               da->e[cnt++] = cell[split[c]];
118454e267fSLisandro Dalcin           }
119454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
120454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
121454e267fSLisandro Dalcin               da->e[cnt++] = cell[c];
122454e267fSLisandro Dalcin           }
123454e267fSLisandro Dalcin         }
124454e267fSLisandro Dalcin       }
125454e267fSLisandro Dalcin     }
126454e267fSLisandro Dalcin   }
127454e267fSLisandro Dalcin   *nel = da->ne;
128454e267fSLisandro Dalcin   *nen = nn;
129454e267fSLisandro Dalcin   *e   = da->e;
130454e267fSLisandro Dalcin   PetscFunctionReturn(0);
131454e267fSLisandro Dalcin }
132454e267fSLisandro Dalcin 
1332dde6fd4SLisandro Dalcin /*@C
1342dde6fd4SLisandro Dalcin       DMDASetElementType - Sets the element type to be returned by DMDAGetElements()
1352dde6fd4SLisandro Dalcin 
1362dde6fd4SLisandro Dalcin     Not Collective
1372dde6fd4SLisandro Dalcin 
1382dde6fd4SLisandro Dalcin    Input Parameter:
1392dde6fd4SLisandro Dalcin .     da - the DMDA object
1402dde6fd4SLisandro Dalcin 
1412dde6fd4SLisandro Dalcin    Output Parameters:
1422dde6fd4SLisandro Dalcin .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
1432dde6fd4SLisandro Dalcin 
1442dde6fd4SLisandro Dalcin    Level: intermediate
1452dde6fd4SLisandro Dalcin 
1462dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
1472dde6fd4SLisandro Dalcin @*/
148454e267fSLisandro Dalcin #undef __FUNCT__
1492dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDASetElementType"
1502dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
1512dde6fd4SLisandro Dalcin {
1522dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
1532dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
1542dde6fd4SLisandro Dalcin 
1552dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1562dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1572dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
1582dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
1592dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
1602dde6fd4SLisandro Dalcin     dd->elementtype = etype;
1612dde6fd4SLisandro Dalcin     dd->ne          = 0;
1622dde6fd4SLisandro Dalcin     dd->e           = PETSC_NULL;
1632dde6fd4SLisandro Dalcin   }
1642dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1652dde6fd4SLisandro Dalcin }
1662dde6fd4SLisandro Dalcin 
1672dde6fd4SLisandro Dalcin /*@C
1682dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
1692dde6fd4SLisandro Dalcin 
1702dde6fd4SLisandro Dalcin     Not Collective
1712dde6fd4SLisandro Dalcin 
1722dde6fd4SLisandro Dalcin    Input Parameter:
1732dde6fd4SLisandro Dalcin .     da - the DMDA object
1742dde6fd4SLisandro Dalcin 
1752dde6fd4SLisandro Dalcin    Output Parameters:
1762dde6fd4SLisandro Dalcin .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
1772dde6fd4SLisandro Dalcin 
1782dde6fd4SLisandro Dalcin    Level: intermediate
1792dde6fd4SLisandro Dalcin 
1802dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
1812dde6fd4SLisandro Dalcin @*/
1822dde6fd4SLisandro Dalcin #undef __FUNCT__
1832dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElementType"
1842dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
1852dde6fd4SLisandro Dalcin {
1862dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
187*5fd66863SKarl Rupp 
1882dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1892dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1902dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
1912dde6fd4SLisandro Dalcin   *etype = dd->elementtype;
1922dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1932dde6fd4SLisandro Dalcin }
1942dde6fd4SLisandro Dalcin 
1952dde6fd4SLisandro Dalcin /*@C
1962dde6fd4SLisandro Dalcin       DMDAGetElements - Gets an array containing the indices (in local coordinates)
1972dde6fd4SLisandro Dalcin                  of all the local elements
1982dde6fd4SLisandro Dalcin 
1992dde6fd4SLisandro Dalcin     Not Collective
2002dde6fd4SLisandro Dalcin 
2012dde6fd4SLisandro Dalcin    Input Parameter:
2022dde6fd4SLisandro Dalcin .     dm - the DM object
2032dde6fd4SLisandro Dalcin 
2042dde6fd4SLisandro Dalcin    Output Parameters:
2052dde6fd4SLisandro Dalcin +     nel - number of local elements
2062dde6fd4SLisandro Dalcin .     nen - number of element nodes
207c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2082dde6fd4SLisandro Dalcin 
2092dde6fd4SLisandro Dalcin    Level: intermediate
2102dde6fd4SLisandro Dalcin 
211c2cd2aa3SJed Brown .seealso: DMDAElementType, DMDASetElementType(), DMDARestoreElements(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
2122dde6fd4SLisandro Dalcin @*/
2132dde6fd4SLisandro Dalcin #undef __FUNCT__
2142dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements"
2152dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
216454e267fSLisandro Dalcin {
217454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
218454e267fSLisandro Dalcin   PetscErrorCode ierr;
219*5fd66863SKarl Rupp 
220454e267fSLisandro Dalcin   PetscFunctionBegin;
221454e267fSLisandro Dalcin   if (da->dim==-1) {
222454e267fSLisandro Dalcin     *nel = 0; *nen = 0; *e = PETSC_NULL;
223454e267fSLisandro Dalcin   } else if (da->dim==1) {
2242dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
225454e267fSLisandro Dalcin   } else if (da->dim==2) {
2262dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
227454e267fSLisandro Dalcin   } else if (da->dim==3) {
2282dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
22930729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",da->dim);
230454e267fSLisandro Dalcin 
231454e267fSLisandro Dalcin   PetscFunctionReturn(0);
232454e267fSLisandro Dalcin }
2332dde6fd4SLisandro Dalcin 
2342dde6fd4SLisandro Dalcin #undef __FUNCT__
2352dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDARestoreElements"
2362dde6fd4SLisandro Dalcin /*@C
2372dde6fd4SLisandro Dalcin       DMDARestoreElements - Returns an array containing the indices (in local coordinates)
2382dde6fd4SLisandro Dalcin                  of all the local elements obtained with DMDAGetElements()
2392dde6fd4SLisandro Dalcin 
2402dde6fd4SLisandro Dalcin     Not Collective
2412dde6fd4SLisandro Dalcin 
2422dde6fd4SLisandro Dalcin    Input Parameter:
2432dde6fd4SLisandro Dalcin +     dm - the DM object
2442dde6fd4SLisandro Dalcin .     nel - number of local elements
2452dde6fd4SLisandro Dalcin .     nen - number of element nodes
246c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2472dde6fd4SLisandro Dalcin 
2482dde6fd4SLisandro Dalcin    Level: intermediate
2492dde6fd4SLisandro Dalcin 
2502dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
2512dde6fd4SLisandro Dalcin @*/
2522dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
2532dde6fd4SLisandro Dalcin {
2542dde6fd4SLisandro Dalcin   PetscFunctionBegin;
2552dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2562dde6fd4SLisandro Dalcin   PetscValidIntPointer(nel,2);
2572dde6fd4SLisandro Dalcin   PetscValidIntPointer(nen,3);
2582dde6fd4SLisandro Dalcin   PetscValidPointer(e,4);
2592dde6fd4SLisandro Dalcin   /* XXX */
2602dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
2612dde6fd4SLisandro Dalcin }
262