xref: /petsc/src/dm/impls/da/dageometry.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
157459e9aSMatthew G Knepley #include <petsc-private/daimpl.h>     /*I  "petscdmda.h"   I*/
257459e9aSMatthew G Knepley 
357459e9aSMatthew G Knepley #undef __FUNCT__
457459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureArray_Private"
557459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar **array)
657459e9aSMatthew G Knepley {
757459e9aSMatthew G Knepley   PetscScalar    *a;
857459e9aSMatthew G Knepley   PetscInt       size = 0, dof, off, d, k, i;
957459e9aSMatthew G Knepley   PetscErrorCode ierr;
1057459e9aSMatthew G Knepley 
1157459e9aSMatthew G Knepley   PetscFunctionBegin;
1257459e9aSMatthew G Knepley   for (i = 0; i < nP; ++i) {
1357459e9aSMatthew G Knepley     ierr  = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
1457459e9aSMatthew G Knepley     size += dof;
1557459e9aSMatthew G Knepley   }
16aa1993deSMatthew G Knepley   ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr);
1757459e9aSMatthew G Knepley   for (i = 0, k = 0; i < nP; ++i) {
1857459e9aSMatthew G Knepley     ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
1957459e9aSMatthew G Knepley     ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
2057459e9aSMatthew G Knepley 
21*8865f1eaSKarl Rupp     for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d];
2257459e9aSMatthew G Knepley   }
2357459e9aSMatthew G Knepley   *array = a;
2457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
2557459e9aSMatthew G Knepley }
2657459e9aSMatthew G Knepley 
2757459e9aSMatthew G Knepley #undef __FUNCT__
2857459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
2957459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
3057459e9aSMatthew G Knepley {
3157459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
3257459e9aSMatthew G Knepley   PetscErrorCode ierr;
3357459e9aSMatthew G Knepley 
3457459e9aSMatthew G Knepley   PetscFunctionBegin;
3557459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
3657459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
3757459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
3857459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
3957459e9aSMatthew G Knepley 
40*8865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k];
4157459e9aSMatthew G Knepley     }
4257459e9aSMatthew G Knepley   } else {
4357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
4457459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4557459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
4657459e9aSMatthew G Knepley 
47*8865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k];
4857459e9aSMatthew G Knepley     }
4957459e9aSMatthew G Knepley   }
5057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
5157459e9aSMatthew G Knepley }
5257459e9aSMatthew G Knepley 
5357459e9aSMatthew G Knepley #undef __FUNCT__
54aa1993deSMatthew G Knepley #define __FUNCT__ "GetPointArray_Private"
55aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
56aa1993deSMatthew G Knepley {
57aa1993deSMatthew G Knepley   PetscErrorCode ierr;
58aa1993deSMatthew G Knepley   PetscInt       *work;
59aa1993deSMatthew G Knepley 
60aa1993deSMatthew G Knepley   PetscFunctionBegin;
61aa1993deSMatthew G Knepley   if (rn) *rn = n;
62aa1993deSMatthew G Knepley   if (rpoints) {
63aa1993deSMatthew G Knepley     ierr     = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
64aa1993deSMatthew G Knepley     ierr     = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
65aa1993deSMatthew G Knepley     *rpoints = work;
66aa1993deSMatthew G Knepley   }
67aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
68aa1993deSMatthew G Knepley }
69aa1993deSMatthew G Knepley 
70aa1993deSMatthew G Knepley #undef __FUNCT__
71aa1993deSMatthew G Knepley #define __FUNCT__ "RestorePointArray_Private"
72aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
73aa1993deSMatthew G Knepley {
74aa1993deSMatthew G Knepley   PetscErrorCode ierr;
75aa1993deSMatthew G Knepley 
76aa1993deSMatthew G Knepley   PetscFunctionBegin;
77aa1993deSMatthew G Knepley   if (rn) *rn = 0;
78*8865f1eaSKarl Rupp   if (rpoints) {
79*8865f1eaSKarl Rupp     ierr = DMRestoreWorkArray(dm,*rn, PETSC_INT, (void*) rpoints);CHKERRQ(ierr);
80*8865f1eaSKarl Rupp   }
81aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
82aa1993deSMatthew G Knepley }
83aa1993deSMatthew G Knepley 
84aa1993deSMatthew G Knepley #undef __FUNCT__
85aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosure"
86aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
8757459e9aSMatthew G Knepley {
8857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
8957459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
90aa1993deSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
91aa1993deSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
92aa1993deSMatthew G Knepley   PetscErrorCode ierr;
93aa1993deSMatthew G Knepley 
94aa1993deSMatthew G Knepley   PetscFunctionBegin;
95aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
96aa1993deSMatthew G Knepley   if (n) PetscValidIntPointer(n,4);
97aa1993deSMatthew G Knepley   if (closure) PetscValidPointer(closure, 5);
98aa1993deSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
99aa1993deSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
100aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
101aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
102aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
103aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
104aa1993deSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
105aa1993deSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
106aa1993deSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
107aa1993deSMatthew G Knepley   yfStart = xfEnd;
108aa1993deSMatthew G Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
109aa1993deSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
110aa1993deSMatthew G Knepley     /* Cell */
111201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
112201c01f8SBarry Smith     else if (dim == 2) {
113aa1993deSMatthew G Knepley       /* 4 faces, 4 vertices
114aa1993deSMatthew G Knepley          Bottom-left vertex follows same order as cells
115aa1993deSMatthew G Knepley          Bottom y-face same order as cells
116aa1993deSMatthew G Knepley          Left x-face follows same order as cells
117aa1993deSMatthew G Knepley          We number the quad:
118aa1993deSMatthew G Knepley 
119aa1993deSMatthew G Knepley            8--3--7
120aa1993deSMatthew G Knepley            |     |
121aa1993deSMatthew G Knepley            4  0  2
122aa1993deSMatthew G Knepley            |     |
123aa1993deSMatthew G Knepley            5--1--6
124aa1993deSMatthew G Knepley       */
125aa1993deSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
126aa1993deSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
127aa1993deSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
128aa1993deSMatthew G Knepley       PetscInt yf = c + yfStart;
129da80777bSKarl Rupp       PetscInt points[9];
130da80777bSKarl Rupp 
131da80777bSKarl Rupp       /* Note: initializer list is not computable at compile time, hence we fill the array manually */
132da80777bSKarl 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;
133aa1993deSMatthew G Knepley 
134aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
135201c01f8SBarry Smith     } else SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
136aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
137aa1993deSMatthew G Knepley     /* Vertex */
138aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
139aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
140aa1993deSMatthew G Knepley     /* X Face */
141201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
142201c01f8SBarry Smith     else if (dim == 2) {
143aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
144aa1993deSMatthew G Knepley       PetscInt f = p - xfStart;
145da80777bSKarl Rupp       PetscInt points[3];
146aa1993deSMatthew G Knepley 
147da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
148aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
149aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
150aa1993deSMatthew G Knepley     } else if (dim == 3) {
151aa1993deSMatthew G Knepley       /* 4 vertices */
152aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
153aa1993deSMatthew G Knepley     }
154aa1993deSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
155aa1993deSMatthew G Knepley     /* Y Face */
156201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
157201c01f8SBarry Smith     else if (dim == 2) {
158aa1993deSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
159aa1993deSMatthew G Knepley       PetscInt f = p - yfStart;
160da80777bSKarl Rupp       PetscInt points[3];
161aa1993deSMatthew G Knepley 
162da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2]= f+1;
163aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
164aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
165aa1993deSMatthew G Knepley     } else if (dim == 3) {
166aa1993deSMatthew G Knepley       /* 4 vertices */
167aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
168aa1993deSMatthew G Knepley     }
169aa1993deSMatthew G Knepley   } else {
170aa1993deSMatthew G Knepley     /* Z Face */
171201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
172201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
173201c01f8SBarry Smith     else if (dim == 3) {
174aa1993deSMatthew G Knepley       /* 4 vertices */
175aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
176aa1993deSMatthew G Knepley     }
177aa1993deSMatthew G Knepley   }
178aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
179aa1993deSMatthew G Knepley }
180aa1993deSMatthew G Knepley 
181aa1993deSMatthew G Knepley #undef __FUNCT__
182aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
183aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
184aa1993deSMatthew G Knepley {
185aa1993deSMatthew G Knepley   PetscErrorCode ierr;
186aa1993deSMatthew G Knepley 
187aa1993deSMatthew G Knepley   PetscFunctionBegin;
188aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
189aa1993deSMatthew G Knepley   PetscValidIntPointer(n,4);
190aa1993deSMatthew G Knepley   PetscValidPointer(closure, 5);
191aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
192aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
193aa1993deSMatthew G Knepley }
194aa1993deSMatthew G Knepley 
195aa1993deSMatthew G Knepley #undef __FUNCT__
196aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
197aa1993deSMatthew G Knepley /* If you did not pass PETSC_NULL for 'values', you must call DMDARestoreClosureScalar() */
198aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
199aa1993deSMatthew G Knepley {
200aa1993deSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
201aa1993deSMatthew G Knepley   PetscInt       dim = da->dim;
20257459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
20395babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
20457459e9aSMatthew G Knepley   PetscErrorCode ierr;
20557459e9aSMatthew G Knepley 
20657459e9aSMatthew G Knepley   PetscFunctionBegin;
20757459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
208aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
20957459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
21057459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
21157459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
21257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
21357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
21457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
21557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
21657459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
21757459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
21857459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
21957459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
22095babc02SJed Brown   zfStart = yfEnd;
22157459e9aSMatthew G Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
22257459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
22357459e9aSMatthew G Knepley     /* Cell */
224201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
225201c01f8SBarry Smith     else if (dim == 2) {
22657459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
22757459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
22857459e9aSMatthew G Knepley          Bottom y-face same order as cells
22957459e9aSMatthew G Knepley          Left x-face follows same order as cells
23057459e9aSMatthew G Knepley          We number the quad:
23157459e9aSMatthew G Knepley 
23257459e9aSMatthew G Knepley            8--3--7
23357459e9aSMatthew G Knepley            |     |
23457459e9aSMatthew G Knepley            4  0  2
23557459e9aSMatthew G Knepley            |     |
23657459e9aSMatthew G Knepley            5--1--6
23757459e9aSMatthew G Knepley       */
23857459e9aSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
23957459e9aSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
24057459e9aSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
24157459e9aSMatthew G Knepley       PetscInt yf = c + yfStart;
242da80777bSKarl Rupp       PetscInt points[9];
24357459e9aSMatthew G Knepley 
244da80777bSKarl 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;
24557459e9aSMatthew G Knepley       ierr      = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr);
24657459e9aSMatthew G Knepley     } else {
24757459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
24857459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
24957459e9aSMatthew G Knepley          Back z-face follows same order as cells
25057459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
25157459e9aSMatthew G Knepley          Left x-face follows same order as cells
25257459e9aSMatthew G Knepley 
25357459e9aSMatthew G Knepley               14-----13
25457459e9aSMatthew G Knepley               /|    /|
25557459e9aSMatthew G Knepley              / | 2 / |
25657459e9aSMatthew G Knepley             / 5|  /  |
25757459e9aSMatthew G Knepley           10-----9  4|
25857459e9aSMatthew G Knepley            |  11-|---12
25957459e9aSMatthew G Knepley            |6 /  |  /
26057459e9aSMatthew G Knepley            | /1 3| /
26157459e9aSMatthew G Knepley            |/    |/
26257459e9aSMatthew G Knepley            7-----8
26357459e9aSMatthew G Knepley       */
26457459e9aSMatthew G Knepley       PetscInt c = p - cStart;
265da80777bSKarl Rupp       PetscInt points[15];
266da80777bSKarl Rupp 
267da80777bSKarl 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;
268da80777bSKarl 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;
269da80777bSKarl Rupp       points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
27057459e9aSMatthew G Knepley 
27157459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
27257459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr);
27357459e9aSMatthew G Knepley     }
27457459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
27557459e9aSMatthew G Knepley     /* Vertex */
27657459e9aSMatthew G Knepley     ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr);
27757459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
27857459e9aSMatthew G Knepley     /* X Face */
279201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
280201c01f8SBarry Smith     else if (dim == 2) {
28157459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
28257459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
283da80777bSKarl Rupp       PetscInt points[3];
28457459e9aSMatthew G Knepley 
285da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
28657459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
28757459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
288*8865f1eaSKarl Rupp     } else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
28957459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
29057459e9aSMatthew G Knepley     /* Y Face */
291201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
292201c01f8SBarry Smith     else if (dim == 2) {
29357459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
29457459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
295da80777bSKarl Rupp       PetscInt points[3];
29657459e9aSMatthew G Knepley 
297da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
29857459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
29957459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
300*8865f1eaSKarl Rupp     } else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
30157459e9aSMatthew G Knepley   } else {
30257459e9aSMatthew G Knepley     /* Z Face */
303201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
304201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
305*8865f1eaSKarl Rupp     else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
30657459e9aSMatthew G Knepley   }
307aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
308aa1993deSMatthew G Knepley }
309aa1993deSMatthew G Knepley 
310aa1993deSMatthew G Knepley #undef __FUNCT__
311aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
312aa1993deSMatthew G Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
313aa1993deSMatthew G Knepley {
314aa1993deSMatthew G Knepley   PetscScalar    *vArray;
315aa1993deSMatthew G Knepley   PetscErrorCode ierr;
316aa1993deSMatthew G Knepley 
317aa1993deSMatthew G Knepley   PetscFunctionBegin;
318aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
319aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
320aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
321aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
322aa1993deSMatthew G Knepley   ierr = DMDAGetClosureScalar(dm,section,p,vArray,values);CHKERRQ(ierr);
32357459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
32457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
32557459e9aSMatthew G Knepley }
32657459e9aSMatthew G Knepley 
32757459e9aSMatthew G Knepley #undef __FUNCT__
328aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
329aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
330aa1993deSMatthew G Knepley {
331aa1993deSMatthew G Knepley   PetscErrorCode ierr;
332aa1993deSMatthew G Knepley   PetscInt       count;
333aa1993deSMatthew G Knepley 
334aa1993deSMatthew G Knepley   PetscFunctionBegin;
335aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
336aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
337aa1993deSMatthew G Knepley   count = 0;                    /* We are lying about the count */
338854432f1SMatthew G Knepley   ierr  = DMRestoreWorkArray(dm, count, PETSC_SCALAR, (void*) values);CHKERRQ(ierr);
339aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
340aa1993deSMatthew G Knepley }
341aa1993deSMatthew G Knepley 
342aa1993deSMatthew G Knepley #undef __FUNCT__
343aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
344aa1993deSMatthew G Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
345aa1993deSMatthew G Knepley {
346aa1993deSMatthew G Knepley   PetscErrorCode ierr;
347aa1993deSMatthew G Knepley 
348aa1993deSMatthew G Knepley   PetscFunctionBegin;
349aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
351aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
352aa1993deSMatthew G Knepley   ierr = DMDARestoreClosureScalar(dm,section,p,PETSC_NULL,values);CHKERRQ(ierr);
353aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
354aa1993deSMatthew G Knepley }
355aa1993deSMatthew G Knepley 
356aa1993deSMatthew G Knepley #undef __FUNCT__
357aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
358aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
35957459e9aSMatthew G Knepley {
36057459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
36157459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
36257459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
36395babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
36457459e9aSMatthew G Knepley   PetscErrorCode ierr;
36557459e9aSMatthew G Knepley 
36657459e9aSMatthew G Knepley   PetscFunctionBegin;
36757459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
368aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
36957459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
37057459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
37157459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
37257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
37357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
37457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
37557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
37657459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
37757459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
37857459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
37957459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
38095babc02SJed Brown   zfStart = yfEnd;
38157459e9aSMatthew G Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
38257459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
38357459e9aSMatthew G Knepley     /* Cell */
384201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
385201c01f8SBarry Smith     else if (dim == 2) {
38657459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
38757459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
38857459e9aSMatthew G Knepley          Bottom y-face same order as cells
38957459e9aSMatthew G Knepley          Left x-face follows same order as cells
39057459e9aSMatthew G Knepley          We number the quad:
39157459e9aSMatthew G Knepley 
39257459e9aSMatthew G Knepley            8--3--7
39357459e9aSMatthew G Knepley            |     |
39457459e9aSMatthew G Knepley            4  0  2
39557459e9aSMatthew G Knepley            |     |
39657459e9aSMatthew G Knepley            5--1--6
39757459e9aSMatthew G Knepley       */
39857459e9aSMatthew G Knepley       PetscInt c = p - cStart;
399da80777bSKarl Rupp       PetscInt points[9];
40057459e9aSMatthew G Knepley 
401da80777bSKarl 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;
402da80777bSKarl Rupp       points[6] = c+vStart+1; points[7] = c+vStart+nVx+1; points[8] = c+vStart+nVx+0;
40357459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
40457459e9aSMatthew G Knepley     } else {
40557459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
40657459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
40757459e9aSMatthew G Knepley          Back z-face follows same order as cells
40857459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
40957459e9aSMatthew G Knepley          Left x-face follows same order as cells
41057459e9aSMatthew G Knepley 
41157459e9aSMatthew G Knepley               14-----13
41257459e9aSMatthew G Knepley               /|    /|
41357459e9aSMatthew G Knepley              / | 2 / |
41457459e9aSMatthew G Knepley             / 5|  /  |
41557459e9aSMatthew G Knepley           10-----9  4|
41657459e9aSMatthew G Knepley            |  11-|---12
41757459e9aSMatthew G Knepley            |6 /  |  /
41857459e9aSMatthew G Knepley            | /1 3| /
41957459e9aSMatthew G Knepley            |/    |/
42057459e9aSMatthew G Knepley            7-----8
42157459e9aSMatthew G Knepley       */
42257459e9aSMatthew G Knepley       PetscInt c = p - cStart;
423da80777bSKarl Rupp       PetscInt points[15];
42457459e9aSMatthew G Knepley 
425da80777bSKarl 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;
426da80777bSKarl 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;
427da80777bSKarl Rupp       points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
42857459e9aSMatthew G Knepley       ierr       = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
42957459e9aSMatthew G Knepley     }
43057459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
43157459e9aSMatthew G Knepley     /* Vertex */
43257459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
43357459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
43457459e9aSMatthew G Knepley     /* X Face */
435201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
436201c01f8SBarry Smith     else if (dim == 2) {
43757459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
43857459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
439da80777bSKarl Rupp       PetscInt points[3];
44057459e9aSMatthew G Knepley 
441da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
44257459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
443*8865f1eaSKarl Rupp     } else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
44457459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
44557459e9aSMatthew G Knepley     /* Y Face */
446201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
447201c01f8SBarry Smith     else if (dim == 2) {
44857459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
44957459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
450da80777bSKarl Rupp       PetscInt points[3];
45157459e9aSMatthew G Knepley 
452da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
45357459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
454*8865f1eaSKarl Rupp     } else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
45557459e9aSMatthew G Knepley   } else {
45657459e9aSMatthew G Knepley     /* Z Face */
457201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
458201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
459*8865f1eaSKarl Rupp     else if (dim == 3) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
46057459e9aSMatthew G Knepley   }
461aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
462aa1993deSMatthew G Knepley }
463aa1993deSMatthew G Knepley 
464aa1993deSMatthew G Knepley #undef __FUNCT__
465aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
466aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
467aa1993deSMatthew G Knepley {
468aa1993deSMatthew G Knepley   PetscScalar    *vArray;
469aa1993deSMatthew G Knepley   PetscErrorCode ierr;
470aa1993deSMatthew G Knepley 
471aa1993deSMatthew G Knepley   PetscFunctionBegin;
472aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
473aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
474aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
475aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
476aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
47757459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
47857459e9aSMatthew G Knepley   PetscFunctionReturn(0);
47957459e9aSMatthew G Knepley }
48057459e9aSMatthew G Knepley 
481ed4e918cSMatthew G Knepley #undef __FUNCT__
482ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
483ed4e918cSMatthew G Knepley /*@
484f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
485341c9bdaSMatthew G Knepley 
486ed4e918cSMatthew G Knepley   Not Collective
487ed4e918cSMatthew G Knepley 
488ed4e918cSMatthew G Knepley   Input Parameter:
489ed4e918cSMatthew G Knepley + da - the distributed array
490ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
491ed4e918cSMatthew G Knepley 
492ed4e918cSMatthew G Knepley   Output Parameter:
493ed4e918cSMatthew G Knepley . cell - the local cell number
494ed4e918cSMatthew G Knepley 
495ed4e918cSMatthew G Knepley   Level: developer
496ed4e918cSMatthew G Knepley 
497ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
498ed4e918cSMatthew G Knepley @*/
499ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
500341c9bdaSMatthew G Knepley {
501341c9bdaSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
502341c9bdaSMatthew G Knepley   const PetscInt dim = da->dim;
5034e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
504ed4e918cSMatthew 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;
505341c9bdaSMatthew G Knepley 
506341c9bdaSMatthew G Knepley   PetscFunctionBegin;
507ed4e918cSMatthew G Knepley   *cell = -1;
508ed4e918cSMatthew G Knepley   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w);
509ed4e918cSMatthew G Knepley   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye);
510ed4e918cSMatthew G Knepley   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze);
511ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
512ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
513341c9bdaSMatthew G Knepley }
514341c9bdaSMatthew G Knepley 
51557459e9aSMatthew G Knepley #undef __FUNCT__
51657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
51715d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
51857459e9aSMatthew G Knepley {
51957459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
52057459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
52157459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
52257459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
52357459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
52457459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
52557459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
52657459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
52757459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
52857459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
52957459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
53057459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
53157459e9aSMatthew G Knepley   PetscReal         invDet;
53257459e9aSMatthew G Knepley   PetscErrorCode    ierr;
53357459e9aSMatthew G Knepley 
53457459e9aSMatthew G Knepley   PetscFunctionBegin;
53515d2e57cSJed Brown #if defined(PETSC_USE_DEBUG)
53615d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
5374c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
53815d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
53915d2e57cSJed Brown #endif
54015d2e57cSJed Brown   J[0]    = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
54115d2e57cSJed Brown   J[2]    = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
54257459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
54357459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
54457459e9aSMatthew G Knepley   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
54557459e9aSMatthew G Knepley   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
54657459e9aSMatthew G Knepley   ierr    = PetscLogFlops(30);CHKERRQ(ierr);
54757459e9aSMatthew G Knepley   PetscFunctionReturn(0);
54857459e9aSMatthew G Knepley }
54957459e9aSMatthew G Knepley 
55057459e9aSMatthew G Knepley #undef __FUNCT__
55157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
55257459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
55357459e9aSMatthew G Knepley {
55457459e9aSMatthew G Knepley   DM                cdm;
55557459e9aSMatthew G Knepley   Vec               coordinates;
55657459e9aSMatthew G Knepley   const PetscScalar *vertices;
55757459e9aSMatthew G Knepley   PetscInt          dim, d, q;
55857459e9aSMatthew G Knepley   PetscErrorCode    ierr;
55957459e9aSMatthew G Knepley 
56057459e9aSMatthew G Knepley   PetscFunctionBegin;
56157459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
56257459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
5636636e97aSMatthew G Knepley   ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
5646636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
56557459e9aSMatthew G Knepley   ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr);
566*8865f1eaSKarl Rupp   for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]);
56757459e9aSMatthew G Knepley   switch (dim) {
56857459e9aSMatthew G Knepley   case 2:
56957459e9aSMatthew G Knepley     for (q = 0; q < quad->numQuadPoints; ++q) {
57057459e9aSMatthew G Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
57157459e9aSMatthew G Knepley     }
57257459e9aSMatthew G Knepley     break;
57357459e9aSMatthew G Knepley   default:
57457459e9aSMatthew G Knepley     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
57557459e9aSMatthew G Knepley   }
57657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
57757459e9aSMatthew G Knepley }
578