xref: /petsc/src/dm/impls/da/dagetelem.c (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
1454e267fSLisandro Dalcin #define PETSCDM_DLL
2454e267fSLisandro Dalcin 
3454e267fSLisandro Dalcin #include "private/daimpl.h"     /*I  "petscdm.h"   I*/
4454e267fSLisandro Dalcin 
5454e267fSLisandro Dalcin #undef __FUNCT__
6454e267fSLisandro Dalcin #define __FUNCT__ "DMGetElements_DA_1D"
7454e267fSLisandro Dalcin static PetscErrorCode DMGetElements_DA_1D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
8454e267fSLisandro Dalcin {
9454e267fSLisandro Dalcin   PetscErrorCode ierr;
10454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
11454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
12454e267fSLisandro Dalcin   PetscInt       cnt=0;
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__
32454e267fSLisandro Dalcin #define __FUNCT__ "DMGetElements_DA_2D"
33454e267fSLisandro Dalcin static PetscErrorCode DMGetElements_DA_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};
42454e267fSLisandro Dalcin   PetscFunctionBegin;
43454e267fSLisandro Dalcin   if (!da->e) {
44454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2; nn=3;}
45454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=4;}
46454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr);
47454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr);
48454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
49454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
50454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1);
51454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
52454e267fSLisandro Dalcin     for (j=ys; j<ye-1; j++) {
53454e267fSLisandro Dalcin       for (i=xs; i<xe-1; i++) {
54454e267fSLisandro Dalcin         cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs);
55454e267fSLisandro Dalcin         cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs);
56454e267fSLisandro Dalcin         cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs);
57454e267fSLisandro Dalcin         cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs);
58454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_P1) {
59454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
60454e267fSLisandro Dalcin             da->e[cnt++] = cell[split[c]];
61454e267fSLisandro Dalcin         }
62454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
63454e267fSLisandro Dalcin           for (c=0; c<ns*nn; c++)
64454e267fSLisandro Dalcin             da->e[cnt++] = cell[c];
65454e267fSLisandro Dalcin         }
66454e267fSLisandro Dalcin       }
67454e267fSLisandro Dalcin     }
68454e267fSLisandro Dalcin   }
69454e267fSLisandro Dalcin   *nel = da->ne;
70454e267fSLisandro Dalcin   *nen = nn;
71454e267fSLisandro Dalcin   *e   = da->e;
72454e267fSLisandro Dalcin   PetscFunctionReturn(0);
73454e267fSLisandro Dalcin }
74454e267fSLisandro Dalcin 
75454e267fSLisandro Dalcin #undef __FUNCT__
76454e267fSLisandro Dalcin #define __FUNCT__ "DMGetElements_DA_3D"
77454e267fSLisandro Dalcin static PetscErrorCode DMGetElements_DA_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
78454e267fSLisandro Dalcin {
79454e267fSLisandro Dalcin   PetscErrorCode ierr;
80454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
81454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
82454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
83454e267fSLisandro Dalcin   PetscInt       k,zs,ze,Zs,Ze;
84454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[8], ns=6, nn=4;
85454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
86454e267fSLisandro Dalcin                                0,1,7,4,
87454e267fSLisandro Dalcin                                1,2,3,7,
88454e267fSLisandro Dalcin                                1,2,7,6,
89454e267fSLisandro Dalcin                                1,4,5,7,
90454e267fSLisandro Dalcin                                1,5,6,7};
91454e267fSLisandro Dalcin   PetscFunctionBegin;
92454e267fSLisandro Dalcin   if (!da->e) {
93454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6; nn=4;}
94454e267fSLisandro Dalcin     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=8;}
95454e267fSLisandro Dalcin     ierr = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
96454e267fSLisandro Dalcin     ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
97454e267fSLisandro Dalcin     xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
98454e267fSLisandro Dalcin     ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
99454e267fSLisandro Dalcin     ze += zs; Ze += Zs; if (zs != Zs) zs -= 1;
100454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
101454e267fSLisandro Dalcin     ierr = PetscMalloc((1 + nn*da->ne)*sizeof(PetscInt),&da->e);CHKERRQ(ierr);
102454e267fSLisandro Dalcin     for (k=zs; k<ze-1; k++) {
103454e267fSLisandro Dalcin       for (j=ys; j<ye-1; j++) {
104454e267fSLisandro Dalcin         for (i=xs; i<xe-1; i++) {
105454e267fSLisandro Dalcin           cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
106454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
107454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
108454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
109454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
110454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
111454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
112454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
113454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
114454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
115454e267fSLisandro Dalcin               da->e[cnt++] = cell[split[c]];
116454e267fSLisandro Dalcin           }
117454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
118454e267fSLisandro Dalcin             for (c=0; c<ns*nn; c++)
119454e267fSLisandro Dalcin               da->e[cnt++] = cell[c];
120454e267fSLisandro Dalcin           }
121454e267fSLisandro Dalcin         }
122454e267fSLisandro Dalcin       }
123454e267fSLisandro Dalcin     }
124454e267fSLisandro Dalcin   }
125454e267fSLisandro Dalcin   *nel = da->ne;
126454e267fSLisandro Dalcin   *nen = nn;
127454e267fSLisandro Dalcin   *e   = da->e;
128454e267fSLisandro Dalcin   PetscFunctionReturn(0);
129454e267fSLisandro Dalcin }
130454e267fSLisandro Dalcin 
131454e267fSLisandro Dalcin #undef __FUNCT__
132454e267fSLisandro Dalcin #define __FUNCT__ "DMGetElements_DA"
133*7087cfbeSBarry Smith PetscErrorCode  DMGetElements_DA(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
134454e267fSLisandro Dalcin {
135454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
136454e267fSLisandro Dalcin   PetscErrorCode ierr;
137454e267fSLisandro Dalcin   PetscFunctionBegin;
138454e267fSLisandro Dalcin   if (da->dim==-1) {
139454e267fSLisandro Dalcin     *nel = 0; *nen = 0; *e = PETSC_NULL;
140454e267fSLisandro Dalcin   } else if (da->dim==1) {
141454e267fSLisandro Dalcin     ierr = DMGetElements_DA_1D(dm,nel,nen,e);CHKERRQ(ierr);
142454e267fSLisandro Dalcin   } else if (da->dim==2) {
143454e267fSLisandro Dalcin     ierr = DMGetElements_DA_2D(dm,nel,nen,e);CHKERRQ(ierr);
144454e267fSLisandro Dalcin   } else if (da->dim==3) {
145454e267fSLisandro Dalcin     ierr = DMGetElements_DA_3D(dm,nel,nen,e);CHKERRQ(ierr);
146454e267fSLisandro Dalcin   } else {
147454e267fSLisandro Dalcin     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",da->dim);
148454e267fSLisandro Dalcin   }
149454e267fSLisandro Dalcin 
150454e267fSLisandro Dalcin   PetscFunctionReturn(0);
151454e267fSLisandro Dalcin }
152