xref: /petsc/src/dm/impls/da/dageometry.c (revision 0a3ada3917fbddb968227e4e816f8d165e7a40ae)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>     /*I  "petscdmda.h"   I*/
257459e9aSMatthew G Knepley 
357459e9aSMatthew G Knepley #undef __FUNCT__
4*0a3ada39SMatthew G. Knepley #define __FUNCT__ "FillClosureArray_Static"
5*0a3ada39SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Static(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, PetscInt *csize, const PetscScalar **array)
657459e9aSMatthew G Knepley {
757459e9aSMatthew G Knepley   PetscScalar    *a;
8*0a3ada39SMatthew G. Knepley   PetscInt       pStart, pEnd, size = 0, dof, off, d, k, i;
957459e9aSMatthew G Knepley   PetscErrorCode ierr;
1057459e9aSMatthew G Knepley 
1157459e9aSMatthew G Knepley   PetscFunctionBegin;
12*0a3ada39SMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
1357459e9aSMatthew G Knepley   for (i = 0; i < nP; ++i) {
14*0a3ada39SMatthew G. Knepley     const PetscInt p = points[i];
15*0a3ada39SMatthew G. Knepley 
16*0a3ada39SMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
17*0a3ada39SMatthew G. Knepley     ierr  = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1857459e9aSMatthew G Knepley     size += dof;
1957459e9aSMatthew G Knepley   }
20*0a3ada39SMatthew G. Knepley   if (csize) *csize = size;
21*0a3ada39SMatthew G. Knepley   if (array) {
22aa1993deSMatthew G Knepley     ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr);
2357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
24*0a3ada39SMatthew G. Knepley       const PetscInt p = points[i];
25*0a3ada39SMatthew G. Knepley 
26*0a3ada39SMatthew G. Knepley       if ((p < pStart) || (p >= pEnd)) continue;
27*0a3ada39SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
28*0a3ada39SMatthew G. Knepley       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
2957459e9aSMatthew G Knepley 
308865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d];
3157459e9aSMatthew G Knepley     }
3257459e9aSMatthew G Knepley     *array = a;
33*0a3ada39SMatthew G. Knepley   }
3457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
3557459e9aSMatthew G Knepley }
3657459e9aSMatthew G Knepley 
3757459e9aSMatthew G Knepley #undef __FUNCT__
3857459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
3957459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
4057459e9aSMatthew G Knepley {
4157459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
4257459e9aSMatthew G Knepley   PetscErrorCode ierr;
4357459e9aSMatthew G Knepley 
4457459e9aSMatthew G Knepley   PetscFunctionBegin;
4557459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
4657459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
4757459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4857459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
4957459e9aSMatthew G Knepley 
508865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k];
5157459e9aSMatthew G Knepley     }
5257459e9aSMatthew G Knepley   } else {
5357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
5457459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
5557459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
5657459e9aSMatthew G Knepley 
578865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k];
5857459e9aSMatthew G Knepley     }
5957459e9aSMatthew G Knepley   }
6057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
6157459e9aSMatthew G Knepley }
6257459e9aSMatthew G Knepley 
6357459e9aSMatthew G Knepley #undef __FUNCT__
64aa1993deSMatthew G Knepley #define __FUNCT__ "GetPointArray_Private"
65aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
66aa1993deSMatthew G Knepley {
67aa1993deSMatthew G Knepley   PetscErrorCode ierr;
68aa1993deSMatthew G Knepley   PetscInt       *work;
69aa1993deSMatthew G Knepley 
70aa1993deSMatthew G Knepley   PetscFunctionBegin;
71aa1993deSMatthew G Knepley   if (rn) *rn = n;
72aa1993deSMatthew G Knepley   if (rpoints) {
73aa1993deSMatthew G Knepley     ierr     = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
74aa1993deSMatthew G Knepley     ierr     = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
75aa1993deSMatthew G Knepley     *rpoints = work;
76aa1993deSMatthew G Knepley   }
77aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
78aa1993deSMatthew G Knepley }
79aa1993deSMatthew G Knepley 
80aa1993deSMatthew G Knepley #undef __FUNCT__
81aa1993deSMatthew G Knepley #define __FUNCT__ "RestorePointArray_Private"
82aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
83aa1993deSMatthew G Knepley {
84aa1993deSMatthew G Knepley   PetscErrorCode ierr;
85aa1993deSMatthew G Knepley 
86aa1993deSMatthew G Knepley   PetscFunctionBegin;
87aa1993deSMatthew G Knepley   if (rn) *rn = 0;
888865f1eaSKarl Rupp   if (rpoints) {
898865f1eaSKarl Rupp     ierr = DMRestoreWorkArray(dm,*rn, PETSC_INT, (void*) rpoints);CHKERRQ(ierr);
908865f1eaSKarl Rupp   }
91aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
92aa1993deSMatthew G Knepley }
93aa1993deSMatthew G Knepley 
94aa1993deSMatthew G Knepley #undef __FUNCT__
95aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosure"
96aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
9757459e9aSMatthew G Knepley {
9857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
9957459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
100aa1993deSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
101aa1993deSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
102aa1993deSMatthew G Knepley   PetscErrorCode ierr;
103aa1993deSMatthew G Knepley 
104aa1993deSMatthew G Knepley   PetscFunctionBegin;
105aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
106aa1993deSMatthew G Knepley   if (n) PetscValidIntPointer(n,4);
107aa1993deSMatthew G Knepley   if (closure) PetscValidPointer(closure, 5);
108aa1993deSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
10982f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
110aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
111aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
112aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
113aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
1140298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
115aa1993deSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
116aa1993deSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
117aa1993deSMatthew G Knepley   yfStart = xfEnd;
11882f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
119aa1993deSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
120aa1993deSMatthew G Knepley     /* Cell */
12182f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
122201c01f8SBarry Smith     else if (dim == 2) {
123aa1993deSMatthew G Knepley       /* 4 faces, 4 vertices
124aa1993deSMatthew G Knepley          Bottom-left vertex follows same order as cells
125aa1993deSMatthew G Knepley          Bottom y-face same order as cells
126aa1993deSMatthew G Knepley          Left x-face follows same order as cells
127aa1993deSMatthew G Knepley          We number the quad:
128aa1993deSMatthew G Knepley 
129aa1993deSMatthew G Knepley            8--3--7
130aa1993deSMatthew G Knepley            |     |
131aa1993deSMatthew G Knepley            4  0  2
132aa1993deSMatthew G Knepley            |     |
133aa1993deSMatthew G Knepley            5--1--6
134aa1993deSMatthew G Knepley       */
135aa1993deSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
136aa1993deSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
137aa1993deSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
138aa1993deSMatthew G Knepley       PetscInt yf = c + yfStart;
139da80777bSKarl Rupp       PetscInt points[9];
140da80777bSKarl Rupp 
141da80777bSKarl Rupp       /* Note: initializer list is not computable at compile time, hence we fill the array manually */
142da80777bSKarl Rupp       points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6]= v+1; points[7] = v+nVx+1; points[8] = v+nVx+0;
143aa1993deSMatthew G Knepley 
144aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
14582f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
146aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
147aa1993deSMatthew G Knepley     /* Vertex */
148aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
149aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
150aa1993deSMatthew G Knepley     /* X Face */
15182f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
152201c01f8SBarry Smith     else if (dim == 2) {
153aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
154aa1993deSMatthew G Knepley       PetscInt f = p - xfStart;
155da80777bSKarl Rupp       PetscInt points[3];
156aa1993deSMatthew G Knepley 
157da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
15882f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
159aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
160aa1993deSMatthew G Knepley     } else if (dim == 3) {
161aa1993deSMatthew G Knepley       /* 4 vertices */
16282f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
163aa1993deSMatthew G Knepley     }
164aa1993deSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
165aa1993deSMatthew G Knepley     /* Y Face */
16682f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
167201c01f8SBarry Smith     else if (dim == 2) {
168aa1993deSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
169aa1993deSMatthew G Knepley       PetscInt f = p - yfStart;
170da80777bSKarl Rupp       PetscInt points[3];
171aa1993deSMatthew G Knepley 
172da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2]= f+1;
17382f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
174aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
175aa1993deSMatthew G Knepley     } else if (dim == 3) {
176aa1993deSMatthew G Knepley       /* 4 vertices */
17782f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
178aa1993deSMatthew G Knepley     }
179aa1993deSMatthew G Knepley   } else {
180aa1993deSMatthew G Knepley     /* Z Face */
18182f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
18282f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
183201c01f8SBarry Smith     else if (dim == 3) {
184aa1993deSMatthew G Knepley       /* 4 vertices */
18582f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
186aa1993deSMatthew G Knepley     }
187aa1993deSMatthew G Knepley   }
188aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
189aa1993deSMatthew G Knepley }
190aa1993deSMatthew G Knepley 
191aa1993deSMatthew G Knepley #undef __FUNCT__
192aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
193aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
194aa1993deSMatthew G Knepley {
195aa1993deSMatthew G Knepley   PetscErrorCode ierr;
196aa1993deSMatthew G Knepley 
197aa1993deSMatthew G Knepley   PetscFunctionBegin;
198aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
199aa1993deSMatthew G Knepley   PetscValidIntPointer(n,4);
200aa1993deSMatthew G Knepley   PetscValidPointer(closure, 5);
201aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
202aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
203aa1993deSMatthew G Knepley }
204aa1993deSMatthew G Knepley 
205aa1993deSMatthew G Knepley #undef __FUNCT__
206aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
207*0a3ada39SMatthew G. Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values)
208aa1993deSMatthew G Knepley {
209aa1993deSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
210aa1993deSMatthew G Knepley   PetscInt       dim = da->dim;
21157459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
21295babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
21357459e9aSMatthew G Knepley   PetscErrorCode ierr;
21457459e9aSMatthew G Knepley 
21557459e9aSMatthew G Knepley   PetscFunctionBegin;
21657459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
217aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
218*0a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
21957459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
22082f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
22157459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
22257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
22357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
22457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
2250298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
22657459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
22757459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
22857459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
22995babc02SJed Brown   zfStart = yfEnd;
23082f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
23157459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
23257459e9aSMatthew G Knepley     /* Cell */
23382f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
234201c01f8SBarry Smith     else if (dim == 2) {
23557459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
23657459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
23757459e9aSMatthew G Knepley          Bottom y-face same order as cells
23857459e9aSMatthew G Knepley          Left x-face follows same order as cells
23957459e9aSMatthew G Knepley          We number the quad:
24057459e9aSMatthew G Knepley 
24157459e9aSMatthew G Knepley            8--3--7
24257459e9aSMatthew G Knepley            |     |
24357459e9aSMatthew G Knepley            4  0  2
24457459e9aSMatthew G Knepley            |     |
24557459e9aSMatthew G Knepley            5--1--6
24657459e9aSMatthew G Knepley       */
24757459e9aSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
24857459e9aSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
24957459e9aSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
25057459e9aSMatthew G Knepley       PetscInt yf = c + yfStart;
251da80777bSKarl Rupp       PetscInt points[9];
25257459e9aSMatthew G Knepley 
253da80777bSKarl Rupp       points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6] = v+1; points[7] = v+nVx+1; points[8] = v+nVx+0;
254*0a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 9, points, vArray, csize, values);CHKERRQ(ierr);
25557459e9aSMatthew G Knepley     } else {
25657459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
25757459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
25857459e9aSMatthew G Knepley          Back z-face follows same order as cells
25957459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
26057459e9aSMatthew G Knepley          Left x-face follows same order as cells
26157459e9aSMatthew G Knepley 
26257459e9aSMatthew G Knepley               14-----13
26357459e9aSMatthew G Knepley               /|    /|
26457459e9aSMatthew G Knepley              / | 2 / |
26557459e9aSMatthew G Knepley             / 5|  /  |
26657459e9aSMatthew G Knepley           10-----9  4|
26757459e9aSMatthew G Knepley            |  11-|---12
26857459e9aSMatthew G Knepley            |6 /  |  /
26957459e9aSMatthew G Knepley            | /1 3| /
27057459e9aSMatthew G Knepley            |/    |/
27157459e9aSMatthew G Knepley            7-----8
27257459e9aSMatthew G Knepley       */
27357459e9aSMatthew G Knepley       PetscInt c = p - cStart;
274da80777bSKarl Rupp       PetscInt points[15];
275da80777bSKarl Rupp 
276da80777bSKarl Rupp       points[0]  = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart;
277da80777bSKarl Rupp       points[7]  = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0;
278da80777bSKarl Rupp       points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
27957459e9aSMatthew G Knepley 
28082f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
281*0a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 15, points, vArray, csize, values);CHKERRQ(ierr);
28257459e9aSMatthew G Knepley     }
28357459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
28457459e9aSMatthew G Knepley     /* Vertex */
285*0a3ada39SMatthew G. Knepley     ierr = FillClosureArray_Static(dm, section, 1, &p, vArray, csize, values);CHKERRQ(ierr);
28657459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
28757459e9aSMatthew G Knepley     /* X Face */
28882f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
289201c01f8SBarry Smith     else if (dim == 2) {
29057459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
29157459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
292da80777bSKarl Rupp       PetscInt points[3];
29357459e9aSMatthew G Knepley 
294da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
29582f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
296*0a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
29782f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
29857459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
29957459e9aSMatthew G Knepley     /* Y Face */
30082f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
301201c01f8SBarry Smith     else if (dim == 2) {
30257459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
30357459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
304da80777bSKarl Rupp       PetscInt points[3];
30557459e9aSMatthew G Knepley 
306da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
30782f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
308*0a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
30982f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
31057459e9aSMatthew G Knepley   } else {
31157459e9aSMatthew G Knepley     /* Z Face */
31282f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
31382f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
31482f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
31557459e9aSMatthew G Knepley   }
316aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
317aa1993deSMatthew G Knepley }
318aa1993deSMatthew G Knepley 
319aa1993deSMatthew G Knepley #undef __FUNCT__
320aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
321*0a3ada39SMatthew G. Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values)
322aa1993deSMatthew G Knepley {
323aa1993deSMatthew G Knepley   PetscScalar    *vArray;
324aa1993deSMatthew G Knepley   PetscErrorCode ierr;
325aa1993deSMatthew G Knepley 
326aa1993deSMatthew G Knepley   PetscFunctionBegin;
327aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
328aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
329*0a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
330aa1993deSMatthew G Knepley   ierr = VecGetArray(v, &vArray);CHKERRQ(ierr);
331*0a3ada39SMatthew G. Knepley   ierr = DMDAGetClosureScalar(dm, section, p, vArray, csize, values);CHKERRQ(ierr);
33257459e9aSMatthew G Knepley   ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr);
33357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
33457459e9aSMatthew G Knepley }
33557459e9aSMatthew G Knepley 
33657459e9aSMatthew G Knepley #undef __FUNCT__
337aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
338*0a3ada39SMatthew G. Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values)
339aa1993deSMatthew G Knepley {
340aa1993deSMatthew G Knepley   PetscErrorCode ierr;
341aa1993deSMatthew G Knepley 
342aa1993deSMatthew G Knepley   PetscFunctionBegin;
343aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
344*0a3ada39SMatthew G. Knepley   PetscValidPointer(values, 6);
345*0a3ada39SMatthew G. Knepley   ierr  = DMRestoreWorkArray(dm, csize ? *csize : 0, PETSC_SCALAR, (void*) values);CHKERRQ(ierr);
346aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
347aa1993deSMatthew G Knepley }
348aa1993deSMatthew G Knepley 
349aa1993deSMatthew G Knepley #undef __FUNCT__
350aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
351*0a3ada39SMatthew G. Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values)
352aa1993deSMatthew G Knepley {
353aa1993deSMatthew G Knepley   PetscErrorCode ierr;
354aa1993deSMatthew G Knepley 
355aa1993deSMatthew G Knepley   PetscFunctionBegin;
356aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
357aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
358aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
359*0a3ada39SMatthew G. Knepley   ierr = DMDARestoreClosureScalar(dm, section, p, NULL, csize, values);CHKERRQ(ierr);
360aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
361aa1993deSMatthew G Knepley }
362aa1993deSMatthew G Knepley 
363aa1993deSMatthew G Knepley #undef __FUNCT__
364aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
365aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
36657459e9aSMatthew G Knepley {
36757459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
36857459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
36957459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
37095babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
37157459e9aSMatthew G Knepley   PetscErrorCode ierr;
37257459e9aSMatthew G Knepley 
37357459e9aSMatthew G Knepley   PetscFunctionBegin;
37457459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
375aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
37657459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
37757459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
37882f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
37957459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
38057459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
38157459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
38257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
3830298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
38457459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
38557459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
38657459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
38795babc02SJed Brown   zfStart = yfEnd;
38882f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
38957459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
39057459e9aSMatthew G Knepley     /* Cell */
39182f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
392201c01f8SBarry Smith     else if (dim == 2) {
39357459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
39457459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
39557459e9aSMatthew G Knepley          Bottom y-face same order as cells
39657459e9aSMatthew G Knepley          Left x-face follows same order as cells
39757459e9aSMatthew G Knepley          We number the quad:
39857459e9aSMatthew G Knepley 
39957459e9aSMatthew G Knepley            8--3--7
40057459e9aSMatthew G Knepley            |     |
40157459e9aSMatthew G Knepley            4  0  2
40257459e9aSMatthew G Knepley            |     |
40357459e9aSMatthew G Knepley            5--1--6
40457459e9aSMatthew G Knepley       */
40557459e9aSMatthew G Knepley       PetscInt c = p - cStart;
406da80777bSKarl Rupp       PetscInt points[9];
40757459e9aSMatthew G Knepley 
408da80777bSKarl Rupp       points[0] = p; points[1] = c+yfStart; points[2] = c+xfStart+1; points[3] = c+yfStart+nyF; points[4] = c+xfStart+0; points[5] = c+vStart+0;
409da80777bSKarl Rupp       points[6] = c+vStart+1; points[7] = c+vStart+nVx+1; points[8] = c+vStart+nVx+0;
41057459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
41157459e9aSMatthew G Knepley     } else {
41257459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
41357459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
41457459e9aSMatthew G Knepley          Back z-face follows same order as cells
41557459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
41657459e9aSMatthew G Knepley          Left x-face follows same order as cells
41757459e9aSMatthew G Knepley 
41857459e9aSMatthew G Knepley               14-----13
41957459e9aSMatthew G Knepley               /|    /|
42057459e9aSMatthew G Knepley              / | 2 / |
42157459e9aSMatthew G Knepley             / 5|  /  |
42257459e9aSMatthew G Knepley           10-----9  4|
42357459e9aSMatthew G Knepley            |  11-|---12
42457459e9aSMatthew G Knepley            |6 /  |  /
42557459e9aSMatthew G Knepley            | /1 3| /
42657459e9aSMatthew G Knepley            |/    |/
42757459e9aSMatthew G Knepley            7-----8
42857459e9aSMatthew G Knepley       */
42957459e9aSMatthew G Knepley       PetscInt c = p - cStart;
430da80777bSKarl Rupp       PetscInt points[15];
43157459e9aSMatthew G Knepley 
432da80777bSKarl Rupp       points[0]  = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart;
433da80777bSKarl Rupp       points[7]  = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; points[12] = c+vStart+nVx*nVy+1;
434da80777bSKarl Rupp       points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
43557459e9aSMatthew G Knepley       ierr       = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
43657459e9aSMatthew G Knepley     }
43757459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
43857459e9aSMatthew G Knepley     /* Vertex */
43957459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
44057459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
44157459e9aSMatthew G Knepley     /* X Face */
44282f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
443201c01f8SBarry Smith     else if (dim == 2) {
44457459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
44557459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
446da80777bSKarl Rupp       PetscInt points[3];
44757459e9aSMatthew G Knepley 
448da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
44957459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
45082f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
45157459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
45257459e9aSMatthew G Knepley     /* Y Face */
45382f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
454201c01f8SBarry Smith     else if (dim == 2) {
45557459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
45657459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
457da80777bSKarl Rupp       PetscInt points[3];
45857459e9aSMatthew G Knepley 
459da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
46057459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
46182f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
46257459e9aSMatthew G Knepley   } else {
46357459e9aSMatthew G Knepley     /* Z Face */
46482f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
46582f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
46682f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
46757459e9aSMatthew G Knepley   }
468aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
469aa1993deSMatthew G Knepley }
470aa1993deSMatthew G Knepley 
471aa1993deSMatthew G Knepley #undef __FUNCT__
472aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
473aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
474aa1993deSMatthew G Knepley {
475aa1993deSMatthew G Knepley   PetscScalar    *vArray;
476aa1993deSMatthew G Knepley   PetscErrorCode ierr;
477aa1993deSMatthew G Knepley 
478aa1993deSMatthew G Knepley   PetscFunctionBegin;
479aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
480aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
481aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
482aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
483aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
48457459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
48557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
48657459e9aSMatthew G Knepley }
48757459e9aSMatthew G Knepley 
488ed4e918cSMatthew G Knepley #undef __FUNCT__
489ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
490ed4e918cSMatthew G Knepley /*@
491f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
492341c9bdaSMatthew G Knepley 
493ed4e918cSMatthew G Knepley   Not Collective
494ed4e918cSMatthew G Knepley 
495ed4e918cSMatthew G Knepley   Input Parameter:
496ed4e918cSMatthew G Knepley + da - the distributed array
497ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
498ed4e918cSMatthew G Knepley 
499ed4e918cSMatthew G Knepley   Output Parameter:
500ed4e918cSMatthew G Knepley . cell - the local cell number
501ed4e918cSMatthew G Knepley 
502ed4e918cSMatthew G Knepley   Level: developer
503ed4e918cSMatthew G Knepley 
504ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
505ed4e918cSMatthew G Knepley @*/
506ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
507341c9bdaSMatthew G Knepley {
508341c9bdaSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
509341c9bdaSMatthew G Knepley   const PetscInt dim = da->dim;
5104e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
511ed4e918cSMatthew G Knepley   const PetscInt il  = s.i - da->Xs/da->w, jl = dim > 1 ? s.j - da->Ys : 0, kl = dim > 2 ? s.k - da->Zs : 0;
512341c9bdaSMatthew G Knepley 
513341c9bdaSMatthew G Knepley   PetscFunctionBegin;
514ed4e918cSMatthew G Knepley   *cell = -1;
51582f516ccSBarry Smith   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w);
51682f516ccSBarry Smith   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye);
51782f516ccSBarry Smith   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze);
518ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
519ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
520341c9bdaSMatthew G Knepley }
521341c9bdaSMatthew G Knepley 
52257459e9aSMatthew G Knepley #undef __FUNCT__
52357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
52415d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
52557459e9aSMatthew G Knepley {
52657459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
52757459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
52857459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
52957459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
53057459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
53157459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
53257459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
53357459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
53457459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
53557459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
53657459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
53757459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
53857459e9aSMatthew G Knepley   PetscReal         invDet;
53957459e9aSMatthew G Knepley   PetscErrorCode    ierr;
54057459e9aSMatthew G Knepley 
54157459e9aSMatthew G Knepley   PetscFunctionBegin;
54215d2e57cSJed Brown #if defined(PETSC_USE_DEBUG)
54315d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
5444c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
54515d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
54615d2e57cSJed Brown #endif
54715d2e57cSJed Brown   J[0]    = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
54815d2e57cSJed Brown   J[2]    = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
54957459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
55057459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
55157459e9aSMatthew G Knepley   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
55257459e9aSMatthew G Knepley   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
55357459e9aSMatthew G Knepley   ierr    = PetscLogFlops(30);CHKERRQ(ierr);
55457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
55557459e9aSMatthew G Knepley }
55657459e9aSMatthew G Knepley 
55757459e9aSMatthew G Knepley #undef __FUNCT__
55857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
55957459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
56057459e9aSMatthew G Knepley {
56157459e9aSMatthew G Knepley   DM                 cdm;
56257459e9aSMatthew G Knepley   Vec                coordinates;
563a1e44745SMatthew G. Knepley   const PetscScalar *vertices = NULL;
564*0a3ada39SMatthew G. Knepley   PetscInt           csize, dim, d, q;
56557459e9aSMatthew G Knepley   PetscErrorCode     ierr;
56657459e9aSMatthew G Knepley 
56757459e9aSMatthew G Knepley   PetscFunctionBegin;
56857459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
56957459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
5706636e97aSMatthew G Knepley   ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
5716636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
572*0a3ada39SMatthew G. Knepley   ierr = DMDAVecGetClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
5738865f1eaSKarl Rupp   for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]);
57457459e9aSMatthew G Knepley   switch (dim) {
57557459e9aSMatthew G Knepley   case 2:
576f9fd7fdbSMatthew G. Knepley     for (q = 0; q < quad->numPoints; ++q) {
577f9fd7fdbSMatthew G. Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->points[q*dim], J, invJ, detJ);CHKERRQ(ierr);
57857459e9aSMatthew G Knepley     }
57957459e9aSMatthew G Knepley     break;
58057459e9aSMatthew G Knepley   default:
58182f516ccSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension %d not supported", dim);
58257459e9aSMatthew G Knepley   }
583*0a3ada39SMatthew G. Knepley   ierr = DMDAVecRestoreClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
58457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
58557459e9aSMatthew G Knepley }
586