xref: /petsc/src/dm/impls/da/dagetelem.c (revision f5f57ec0aaa1e5be6a6f3c2ef714a1a38ff0b64a)
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) {
1398b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
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);
18854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + 2*da->ne,&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 
302dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
31454e267fSLisandro Dalcin {
32454e267fSLisandro Dalcin   PetscErrorCode ierr;
33454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
34454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
35454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
36454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[4], ns=2, nn=3;
37454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,
38454e267fSLisandro Dalcin                                2,3,1};
395fd66863SKarl Rupp 
40454e267fSLisandro Dalcin   PetscFunctionBegin;
41b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_P1) {nn=3;}
42b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_Q1) {nn=4;}
43454e267fSLisandro Dalcin   if (!da->e) {
4498b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
45b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=2;}
46b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
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);
52854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&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) {
608865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]];
61454e267fSLisandro Dalcin         }
62454e267fSLisandro Dalcin         if (da->elementtype == DMDA_ELEMENT_Q1) {
638865f1eaSKarl Rupp           for (c=0; c<ns*nn; c++) 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 
742dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
75454e267fSLisandro Dalcin {
76454e267fSLisandro Dalcin   PetscErrorCode ierr;
77454e267fSLisandro Dalcin   DM_DA          *da = (DM_DA*)dm->data;
78454e267fSLisandro Dalcin   PetscInt       i,xs,xe,Xs,Xe;
79454e267fSLisandro Dalcin   PetscInt       j,ys,ye,Ys,Ye;
80454e267fSLisandro Dalcin   PetscInt       k,zs,ze,Zs,Ze;
81454e267fSLisandro Dalcin   PetscInt       cnt=0, cell[8], ns=6, nn=4;
82454e267fSLisandro Dalcin   PetscInt       c, split[] = {0,1,3,7,
83454e267fSLisandro Dalcin                                0,1,7,4,
84454e267fSLisandro Dalcin                                1,2,3,7,
85454e267fSLisandro Dalcin                                1,2,7,6,
86454e267fSLisandro Dalcin                                1,4,5,7,
87454e267fSLisandro Dalcin                                1,5,6,7};
885fd66863SKarl Rupp 
89454e267fSLisandro Dalcin   PetscFunctionBegin;
90b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_P1) {nn=4;}
91b644c393SDave May   if (da->elementtype == DMDA_ELEMENT_Q1) {nn=8;}
92454e267fSLisandro Dalcin   if (!da->e) {
9398b2d131SBarry Smith     if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width");
94b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_P1) {ns=6;}
95b644c393SDave May     if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;}
96454e267fSLisandro Dalcin     ierr   = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
97454e267fSLisandro Dalcin     ierr   = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
98454e267fSLisandro Dalcin     xe    += xs; Xe += Xs; if (xs != Xs) xs -= 1;
99454e267fSLisandro Dalcin     ye    += ys; Ye += Ys; if (ys != Ys) ys -= 1;
100454e267fSLisandro Dalcin     ze    += zs; Ze += Zs; if (zs != Zs) zs -= 1;
101454e267fSLisandro Dalcin     da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1);
102854ce69bSBarry Smith     ierr   = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr);
103454e267fSLisandro Dalcin     for (k=zs; k<ze-1; k++) {
104454e267fSLisandro Dalcin       for (j=ys; j<ye-1; j++) {
105454e267fSLisandro Dalcin         for (i=xs; i<xe-1; i++) {
106454e267fSLisandro Dalcin           cell[0] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
107454e267fSLisandro Dalcin           cell[1] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
108454e267fSLisandro Dalcin           cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
109454e267fSLisandro Dalcin           cell[3] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs  )*(Xe-Xs)*(Ye-Ys);
110454e267fSLisandro Dalcin           cell[4] = (i-Xs  ) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
111454e267fSLisandro Dalcin           cell[5] = (i-Xs+1) + (j-Ys  )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
112454e267fSLisandro Dalcin           cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
113454e267fSLisandro Dalcin           cell[7] = (i-Xs  ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys);
114454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_P1) {
1158865f1eaSKarl Rupp             for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]];
116454e267fSLisandro Dalcin           }
117454e267fSLisandro Dalcin           if (da->elementtype == DMDA_ELEMENT_Q1) {
1188865f1eaSKarl Rupp             for (c=0; c<ns*nn; c++) 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 
130*f5f57ec0SBarry Smith /*@
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:
1396420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
1402dde6fd4SLisandro Dalcin 
1412dde6fd4SLisandro Dalcin    Level: intermediate
1422dde6fd4SLisandro Dalcin 
1432dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements()
1442dde6fd4SLisandro Dalcin @*/
1452dde6fd4SLisandro Dalcin PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
1462dde6fd4SLisandro Dalcin {
1472dde6fd4SLisandro Dalcin   DM_DA          *dd = (DM_DA*)da->data;
1482dde6fd4SLisandro Dalcin   PetscErrorCode ierr;
1492dde6fd4SLisandro Dalcin 
1502dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1512dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1522dde6fd4SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,etype,2);
1532dde6fd4SLisandro Dalcin   if (dd->elementtype != etype) {
1542dde6fd4SLisandro Dalcin     ierr = PetscFree(dd->e);CHKERRQ(ierr);
1558865f1eaSKarl Rupp 
1562dde6fd4SLisandro Dalcin     dd->elementtype = etype;
1572dde6fd4SLisandro Dalcin     dd->ne          = 0;
1580298fd71SBarry Smith     dd->e           = NULL;
1592dde6fd4SLisandro Dalcin   }
1602dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1612dde6fd4SLisandro Dalcin }
1622dde6fd4SLisandro Dalcin 
163*f5f57ec0SBarry Smith /*@
1642dde6fd4SLisandro Dalcin       DMDAGetElementType - Gets the element type to be returned by DMDAGetElements()
1652dde6fd4SLisandro Dalcin 
1662dde6fd4SLisandro Dalcin     Not Collective
1672dde6fd4SLisandro Dalcin 
1682dde6fd4SLisandro Dalcin    Input Parameter:
1692dde6fd4SLisandro Dalcin .     da - the DMDA object
1702dde6fd4SLisandro Dalcin 
1712dde6fd4SLisandro Dalcin    Output Parameters:
1726420f107SDave May .     etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1
1732dde6fd4SLisandro Dalcin 
1742dde6fd4SLisandro Dalcin    Level: intermediate
1752dde6fd4SLisandro Dalcin 
1762dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements()
1772dde6fd4SLisandro Dalcin @*/
1782dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
1792dde6fd4SLisandro Dalcin {
1802dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
1815fd66863SKarl Rupp 
1822dde6fd4SLisandro Dalcin   PetscFunctionBegin;
1832dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1842dde6fd4SLisandro Dalcin   PetscValidPointer(etype,2);
1852dde6fd4SLisandro Dalcin   *etype = dd->elementtype;
1862dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
1872dde6fd4SLisandro Dalcin }
1882dde6fd4SLisandro Dalcin 
1892dde6fd4SLisandro Dalcin /*@C
1902dde6fd4SLisandro Dalcin       DMDAGetElements - Gets an array containing the indices (in local coordinates)
1912dde6fd4SLisandro Dalcin                  of all the local elements
1922dde6fd4SLisandro Dalcin 
1932dde6fd4SLisandro Dalcin     Not Collective
1942dde6fd4SLisandro Dalcin 
1952dde6fd4SLisandro Dalcin    Input Parameter:
1962dde6fd4SLisandro Dalcin .     dm - the DM object
1972dde6fd4SLisandro Dalcin 
1982dde6fd4SLisandro Dalcin    Output Parameters:
1992dde6fd4SLisandro Dalcin +     nel - number of local elements
2002dde6fd4SLisandro Dalcin .     nen - number of element nodes
201c2cd2aa3SJed Brown -     e - the local indices of the elements' vertices
2022dde6fd4SLisandro Dalcin 
2032dde6fd4SLisandro Dalcin    Level: intermediate
2042dde6fd4SLisandro Dalcin 
20593386343SJed Brown    Notes:
20693386343SJed Brown      Call DMDARestoreElements() once you have finished accessing the elements.
20793386343SJed Brown 
20845950774SSatish Balay      Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes.
209009f0576SBarry Smith 
210009f0576SBarry 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.
211009f0576SBarry Smith 
212*f5f57ec0SBarry Smith      Not supported in Fortran
213*f5f57ec0SBarry Smith 
21493386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin()
2152dde6fd4SLisandro Dalcin @*/
2162dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
217454e267fSLisandro Dalcin {
218c73cfb54SMatthew G. Knepley   PetscInt       dim;
219454e267fSLisandro Dalcin   PetscErrorCode ierr;
22099c57625SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
2215fd66863SKarl Rupp 
222454e267fSLisandro Dalcin   PetscFunctionBegin;
2239192b5d6SBarry Smith   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");
224c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
225c73cfb54SMatthew G. Knepley   if (dim==-1) {
2260298fd71SBarry Smith     *nel = 0; *nen = 0; *e = NULL;
227c73cfb54SMatthew G. Knepley   } else if (dim==1) {
2282dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr);
229c73cfb54SMatthew G. Knepley   } else if (dim==2) {
2302dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr);
231c73cfb54SMatthew G. Knepley   } else if (dim==3) {
2322dde6fd4SLisandro Dalcin     ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr);
233c73cfb54SMatthew G. Knepley   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
234454e267fSLisandro Dalcin   PetscFunctionReturn(0);
235454e267fSLisandro Dalcin }
2362dde6fd4SLisandro Dalcin 
2372dde6fd4SLisandro Dalcin /*@C
238009f0576SBarry Smith       DMDARestoreElements - Restores the array 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 
250009f0576SBarry Smith    Note: You should not access these values after you have called this routine.
251009f0576SBarry Smith 
252c67aec63SBarry Smith          This restore signals the DMDA object that you no longer need access to the array information.
253009f0576SBarry Smith 
254*f5f57ec0SBarry Smith          Not supported in Fortran
255*f5f57ec0SBarry Smith 
2562dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements()
2572dde6fd4SLisandro Dalcin @*/
2582dde6fd4SLisandro Dalcin PetscErrorCode  DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
2592dde6fd4SLisandro Dalcin {
2602dde6fd4SLisandro Dalcin   PetscFunctionBegin;
2612dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2622dde6fd4SLisandro Dalcin   PetscValidIntPointer(nel,2);
2632dde6fd4SLisandro Dalcin   PetscValidIntPointer(nen,3);
2642dde6fd4SLisandro Dalcin   PetscValidPointer(e,4);
2659cc8bc54SJed Brown   *nel = 0;
2669cc8bc54SJed Brown   *nen = -1;
2679cc8bc54SJed Brown   *e = NULL;
2682dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
2692dde6fd4SLisandro Dalcin }
270