xref: /petsc/src/dm/impls/da/dageometry.c (revision 201c01f886d7be7e3807209315967bf58c11f9dc)
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 
2157459e9aSMatthew G Knepley     for (d = 0; d < dof; ++d, ++k) {
2257459e9aSMatthew G Knepley       a[k] = vArray[off+d];
2357459e9aSMatthew G Knepley     }
2457459e9aSMatthew G Knepley   }
2557459e9aSMatthew G Knepley   *array = a;
2657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
2757459e9aSMatthew G Knepley }
2857459e9aSMatthew G Knepley 
2957459e9aSMatthew G Knepley #undef __FUNCT__
3057459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
3157459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
3257459e9aSMatthew G Knepley {
3357459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
3457459e9aSMatthew G Knepley   PetscErrorCode ierr;
3557459e9aSMatthew G Knepley 
3657459e9aSMatthew G Knepley   PetscFunctionBegin;
3757459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
3857459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
3957459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4057459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
4157459e9aSMatthew G Knepley 
4257459e9aSMatthew G Knepley       for (d = 0; d < dof; ++d, ++k) {
4357459e9aSMatthew G Knepley         vArray[off+d] = array[k];
4457459e9aSMatthew G Knepley       }
4557459e9aSMatthew G Knepley     }
4657459e9aSMatthew G Knepley   } else {
4757459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
4857459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4957459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
5057459e9aSMatthew G Knepley 
5157459e9aSMatthew G Knepley       for (d = 0; d < dof; ++d, ++k) {
5257459e9aSMatthew G Knepley         vArray[off+d] += array[k];
5357459e9aSMatthew G Knepley       }
5457459e9aSMatthew G Knepley     }
5557459e9aSMatthew G Knepley   }
5657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
5757459e9aSMatthew G Knepley }
5857459e9aSMatthew G Knepley 
5957459e9aSMatthew G Knepley #undef __FUNCT__
60aa1993deSMatthew G Knepley #define __FUNCT__ "GetPointArray_Private"
61aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
62aa1993deSMatthew G Knepley {
63aa1993deSMatthew G Knepley   PetscErrorCode ierr;
64aa1993deSMatthew G Knepley   PetscInt *work;
65aa1993deSMatthew G Knepley 
66aa1993deSMatthew G Knepley   PetscFunctionBegin;
67aa1993deSMatthew G Knepley   if (rn) *rn = n;
68aa1993deSMatthew G Knepley   if (rpoints) {
69aa1993deSMatthew G Knepley     ierr = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
70aa1993deSMatthew G Knepley     ierr = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
71aa1993deSMatthew G Knepley     *rpoints = work;
72aa1993deSMatthew G Knepley   }
73aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
74aa1993deSMatthew G Knepley }
75aa1993deSMatthew G Knepley 
76aa1993deSMatthew G Knepley #undef __FUNCT__
77aa1993deSMatthew G Knepley #define __FUNCT__ "RestorePointArray_Private"
78aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
79aa1993deSMatthew G Knepley {
80aa1993deSMatthew G Knepley   PetscErrorCode ierr;
81aa1993deSMatthew G Knepley 
82aa1993deSMatthew G Knepley   PetscFunctionBegin;
83aa1993deSMatthew G Knepley   if (rn) *rn = 0;
84854432f1SMatthew G Knepley   if (rpoints) {ierr = DMRestoreWorkArray(dm,* rn, PETSC_INT, (void *) rpoints);CHKERRQ(ierr);}
85aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
86aa1993deSMatthew G Knepley }
87aa1993deSMatthew G Knepley 
88aa1993deSMatthew G Knepley #undef __FUNCT__
89aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosure"
90aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
9157459e9aSMatthew G Knepley {
9257459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
9357459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
94aa1993deSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
95aa1993deSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
96aa1993deSMatthew G Knepley   PetscErrorCode ierr;
97aa1993deSMatthew G Knepley 
98aa1993deSMatthew G Knepley   PetscFunctionBegin;
99aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
100aa1993deSMatthew G Knepley   if (n) PetscValidIntPointer(n,4);
101aa1993deSMatthew G Knepley   if (closure) PetscValidPointer(closure, 5);
102aa1993deSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
103aa1993deSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
104aa1993deSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
105aa1993deSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
106aa1993deSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
107aa1993deSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
108aa1993deSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
109aa1993deSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
110aa1993deSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
111aa1993deSMatthew G Knepley   yfStart = xfEnd;
112aa1993deSMatthew 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);
113aa1993deSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
114aa1993deSMatthew G Knepley     /* Cell */
115*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
116*201c01f8SBarry Smith     else if (dim == 2) {
117aa1993deSMatthew G Knepley       /* 4 faces, 4 vertices
118aa1993deSMatthew G Knepley          Bottom-left vertex follows same order as cells
119aa1993deSMatthew G Knepley          Bottom y-face same order as cells
120aa1993deSMatthew G Knepley          Left x-face follows same order as cells
121aa1993deSMatthew G Knepley          We number the quad:
122aa1993deSMatthew G Knepley 
123aa1993deSMatthew G Knepley            8--3--7
124aa1993deSMatthew G Knepley            |     |
125aa1993deSMatthew G Knepley            4  0  2
126aa1993deSMatthew G Knepley            |     |
127aa1993deSMatthew G Knepley            5--1--6
128aa1993deSMatthew G Knepley       */
129aa1993deSMatthew G Knepley       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
130aa1993deSMatthew G Knepley       PetscInt v         = cy*nVx + cx +  vStart;
131aa1993deSMatthew G Knepley       PetscInt xf        = cy*nxF + cx + xfStart;
132aa1993deSMatthew G Knepley       PetscInt yf        = c + yfStart;
133aa1993deSMatthew G Knepley       PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0};
134aa1993deSMatthew G Knepley 
135aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
136*201c01f8SBarry Smith     } else SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
137aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
138aa1993deSMatthew G Knepley     /* Vertex */
139aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
140aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
141aa1993deSMatthew G Knepley     /* X Face */
142*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
143*201c01f8SBarry Smith     else if (dim == 2) {
144aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
145aa1993deSMatthew G Knepley       PetscInt f         = p - xfStart;
146aa1993deSMatthew G Knepley       PetscInt points[3] = {p, f, f+nVx};
147aa1993deSMatthew G Knepley 
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 */
156*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
157*201c01f8SBarry 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;
160aa1993deSMatthew G Knepley       PetscInt points[3] = {p, f, f+1};
161aa1993deSMatthew G Knepley 
162aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
163aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
164aa1993deSMatthew G Knepley     } else if (dim == 3) {
165aa1993deSMatthew G Knepley       /* 4 vertices */
166aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
167aa1993deSMatthew G Knepley     }
168aa1993deSMatthew G Knepley   } else {
169aa1993deSMatthew G Knepley     /* Z Face */
170*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
171*201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
172*201c01f8SBarry Smith     else if (dim == 3) {
173aa1993deSMatthew G Knepley       /* 4 vertices */
174aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
175aa1993deSMatthew G Knepley     }
176aa1993deSMatthew G Knepley   }
177aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
178aa1993deSMatthew G Knepley }
179aa1993deSMatthew G Knepley 
180aa1993deSMatthew G Knepley #undef __FUNCT__
181aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
182aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
183aa1993deSMatthew G Knepley {
184aa1993deSMatthew G Knepley   PetscErrorCode ierr;
185aa1993deSMatthew G Knepley 
186aa1993deSMatthew G Knepley   PetscFunctionBegin;
187aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
188aa1993deSMatthew G Knepley   PetscValidIntPointer(n,4);
189aa1993deSMatthew G Knepley   PetscValidPointer(closure, 5);
190aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
191aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
192aa1993deSMatthew G Knepley }
193aa1993deSMatthew G Knepley 
194aa1993deSMatthew G Knepley #undef __FUNCT__
195aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
196aa1993deSMatthew G Knepley /* If you did not pass PETSC_NULL for 'values', you must call DMDARestoreClosureScalar() */
197aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
198aa1993deSMatthew G Knepley {
199aa1993deSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
200aa1993deSMatthew G Knepley   PetscInt       dim = da->dim;
20157459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
20295babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
20357459e9aSMatthew G Knepley   PetscErrorCode ierr;
20457459e9aSMatthew G Knepley 
20557459e9aSMatthew G Knepley   PetscFunctionBegin;
20657459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
207aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
20857459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
20957459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
21057459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
21157459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
21257459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
21357459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
21457459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
21557459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
21657459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
21757459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
21857459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
21995babc02SJed Brown   zfStart = yfEnd;
22057459e9aSMatthew 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);
22157459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
22257459e9aSMatthew G Knepley     /* Cell */
223*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
224*201c01f8SBarry Smith     else if (dim == 2) {
22557459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
22657459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
22757459e9aSMatthew G Knepley          Bottom y-face same order as cells
22857459e9aSMatthew G Knepley          Left x-face follows same order as cells
22957459e9aSMatthew G Knepley          We number the quad:
23057459e9aSMatthew G Knepley 
23157459e9aSMatthew G Knepley            8--3--7
23257459e9aSMatthew G Knepley            |     |
23357459e9aSMatthew G Knepley            4  0  2
23457459e9aSMatthew G Knepley            |     |
23557459e9aSMatthew G Knepley            5--1--6
23657459e9aSMatthew G Knepley       */
23757459e9aSMatthew G Knepley       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
23857459e9aSMatthew G Knepley       PetscInt v         = cy*nVx + cx +  vStart;
23957459e9aSMatthew G Knepley       PetscInt xf        = cy*nxF + cx + xfStart;
24057459e9aSMatthew G Knepley       PetscInt yf        = c + yfStart;
24157459e9aSMatthew G Knepley       PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0};
24257459e9aSMatthew G Knepley 
24357459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr);
24457459e9aSMatthew G Knepley     } else {
24557459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
24657459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
24757459e9aSMatthew G Knepley          Back z-face follows same order as cells
24857459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
24957459e9aSMatthew G Knepley          Left x-face follows same order as cells
25057459e9aSMatthew G Knepley 
25157459e9aSMatthew G Knepley               14-----13
25257459e9aSMatthew G Knepley               /|    /|
25357459e9aSMatthew G Knepley              / | 2 / |
25457459e9aSMatthew G Knepley             / 5|  /  |
25557459e9aSMatthew G Knepley           10-----9  4|
25657459e9aSMatthew G Knepley            |  11-|---12
25757459e9aSMatthew G Knepley            |6 /  |  /
25857459e9aSMatthew G Knepley            | /1 3| /
25957459e9aSMatthew G Knepley            |/    |/
26057459e9aSMatthew G Knepley            7-----8
26157459e9aSMatthew G Knepley       */
26257459e9aSMatthew G Knepley       PetscInt c          = p - cStart;
26357459e9aSMatthew G Knepley       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
26457459e9aSMatthew G Knepley                              c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0};
26557459e9aSMatthew G Knepley 
26657459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
26757459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr);
26857459e9aSMatthew G Knepley     }
26957459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
27057459e9aSMatthew G Knepley     /* Vertex */
27157459e9aSMatthew G Knepley     ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr);
27257459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
27357459e9aSMatthew G Knepley     /* X Face */
274*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
275*201c01f8SBarry Smith     else if (dim == 2) {
27657459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
27757459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
27857459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+nVx};
27957459e9aSMatthew G Knepley 
28057459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
28157459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
28257459e9aSMatthew G Knepley     } else if (dim == 3) {
28357459e9aSMatthew G Knepley       /* 4 vertices */
28457459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
28557459e9aSMatthew G Knepley     }
28657459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
28757459e9aSMatthew G Knepley     /* Y Face */
288*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
289*201c01f8SBarry Smith     else if (dim == 2) {
29057459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
29157459e9aSMatthew G Knepley       PetscInt f         = p - yfStart;
29257459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+1};
29357459e9aSMatthew G Knepley 
29457459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
29557459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
29657459e9aSMatthew G Knepley     } else if (dim == 3) {
29757459e9aSMatthew G Knepley       /* 4 vertices */
29857459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
29957459e9aSMatthew G Knepley     }
30057459e9aSMatthew G Knepley   } else {
30157459e9aSMatthew G Knepley     /* Z Face */
302*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
303*201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
304*201c01f8SBarry Smith     else if (dim == 3) {
30557459e9aSMatthew G Knepley       /* 4 vertices */
30657459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
30757459e9aSMatthew G Knepley     }
30857459e9aSMatthew G Knepley   }
309aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
310aa1993deSMatthew G Knepley }
311aa1993deSMatthew G Knepley 
312aa1993deSMatthew G Knepley #undef __FUNCT__
313aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
314aa1993deSMatthew G Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
315aa1993deSMatthew G Knepley {
316aa1993deSMatthew G Knepley   PetscScalar *vArray;
317aa1993deSMatthew G Knepley   PetscErrorCode ierr;
318aa1993deSMatthew G Knepley 
319aa1993deSMatthew G Knepley   PetscFunctionBegin;
320aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
321aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
322aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
323aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
324aa1993deSMatthew G Knepley   ierr = DMDAGetClosureScalar(dm,section,p,vArray,values);CHKERRQ(ierr);
32557459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
32657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
32757459e9aSMatthew G Knepley }
32857459e9aSMatthew G Knepley 
32957459e9aSMatthew G Knepley #undef __FUNCT__
330aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
331aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
332aa1993deSMatthew G Knepley {
333aa1993deSMatthew G Knepley   PetscErrorCode ierr;
334aa1993deSMatthew G Knepley   PetscInt count;
335aa1993deSMatthew G Knepley 
336aa1993deSMatthew G Knepley   PetscFunctionBegin;
337aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
338aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
339aa1993deSMatthew G Knepley   count = 0;                    /* We are lying about the count */
340854432f1SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, count, PETSC_SCALAR, (void *) values);CHKERRQ(ierr);
341aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
342aa1993deSMatthew G Knepley }
343aa1993deSMatthew G Knepley 
344aa1993deSMatthew G Knepley #undef __FUNCT__
345aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
346aa1993deSMatthew G Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
347aa1993deSMatthew G Knepley {
348aa1993deSMatthew G Knepley   PetscErrorCode ierr;
349aa1993deSMatthew G Knepley 
350aa1993deSMatthew G Knepley   PetscFunctionBegin;
351aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
352aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
353aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
354aa1993deSMatthew G Knepley   ierr = DMDARestoreClosureScalar(dm,section,p,PETSC_NULL,values);CHKERRQ(ierr);
355aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
356aa1993deSMatthew G Knepley }
357aa1993deSMatthew G Knepley 
358aa1993deSMatthew G Knepley #undef __FUNCT__
359aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
360aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
36157459e9aSMatthew G Knepley {
36257459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
36357459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
36457459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
36595babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
36657459e9aSMatthew G Knepley   PetscErrorCode ierr;
36757459e9aSMatthew G Knepley 
36857459e9aSMatthew G Knepley   PetscFunctionBegin;
36957459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
370aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
37157459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
37257459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
37357459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
37457459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
37557459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
37657459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
37757459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
37857459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
37957459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
38057459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
38157459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
38295babc02SJed Brown   zfStart = yfEnd;
38357459e9aSMatthew 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);
38457459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
38557459e9aSMatthew G Knepley     /* Cell */
386*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
387*201c01f8SBarry Smith     else if (dim == 2) {
38857459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
38957459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
39057459e9aSMatthew G Knepley          Bottom y-face same order as cells
39157459e9aSMatthew G Knepley          Left x-face follows same order as cells
39257459e9aSMatthew G Knepley          We number the quad:
39357459e9aSMatthew G Knepley 
39457459e9aSMatthew G Knepley            8--3--7
39557459e9aSMatthew G Knepley            |     |
39657459e9aSMatthew G Knepley            4  0  2
39757459e9aSMatthew G Knepley            |     |
39857459e9aSMatthew G Knepley            5--1--6
39957459e9aSMatthew G Knepley       */
40057459e9aSMatthew G Knepley       PetscInt c         = p - cStart;
40157459e9aSMatthew G Knepley       PetscInt points[9] = {p, c+yfStart, c+xfStart+1, c+yfStart+nyF, c+xfStart+0, c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0};
40257459e9aSMatthew G Knepley 
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;
42357459e9aSMatthew G Knepley       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
42457459e9aSMatthew G Knepley                              c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0};
42557459e9aSMatthew G Knepley 
42657459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
42757459e9aSMatthew G Knepley     }
42857459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
42957459e9aSMatthew G Knepley     /* Vertex */
43057459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
43157459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
43257459e9aSMatthew G Knepley     /* X Face */
433*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
434*201c01f8SBarry Smith     else if (dim == 2) {
43557459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
43657459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
43757459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+nVx};
43857459e9aSMatthew G Knepley 
43957459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
44057459e9aSMatthew G Knepley     } else if (dim == 3) {
44157459e9aSMatthew G Knepley       /* 4 vertices */
44257459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
44357459e9aSMatthew G Knepley     }
44457459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
44557459e9aSMatthew G Knepley     /* Y Face */
446*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
447*201c01f8SBarry 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;
45057459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+1};
45157459e9aSMatthew G Knepley 
45257459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
45357459e9aSMatthew G Knepley     } else if (dim == 3) {
45457459e9aSMatthew G Knepley       /* 4 vertices */
45557459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
45657459e9aSMatthew G Knepley     }
45757459e9aSMatthew G Knepley   } else {
45857459e9aSMatthew G Knepley     /* Z Face */
459*201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
460*201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
461*201c01f8SBarry Smith     else if (dim == 3) {
46257459e9aSMatthew G Knepley       /* 4 vertices */
46357459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
46457459e9aSMatthew G Knepley     }
46557459e9aSMatthew G Knepley   }
466aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
467aa1993deSMatthew G Knepley }
468aa1993deSMatthew G Knepley 
469aa1993deSMatthew G Knepley #undef __FUNCT__
470aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
471aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
472aa1993deSMatthew G Knepley {
473aa1993deSMatthew G Knepley   PetscScalar   *vArray;
474aa1993deSMatthew G Knepley   PetscErrorCode ierr;
475aa1993deSMatthew G Knepley 
476aa1993deSMatthew G Knepley   PetscFunctionBegin;
477aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
478aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
479aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
480aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
481aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
48257459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
48357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
48457459e9aSMatthew G Knepley }
48557459e9aSMatthew G Knepley 
486ed4e918cSMatthew G Knepley #undef __FUNCT__
487ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
488ed4e918cSMatthew G Knepley /*@
489f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
490341c9bdaSMatthew G Knepley 
491ed4e918cSMatthew G Knepley   Not Collective
492ed4e918cSMatthew G Knepley 
493ed4e918cSMatthew G Knepley   Input Parameter:
494ed4e918cSMatthew G Knepley + da - the distributed array
495ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
496ed4e918cSMatthew G Knepley 
497ed4e918cSMatthew G Knepley   Output Parameter:
498ed4e918cSMatthew G Knepley . cell - the local cell number
499ed4e918cSMatthew G Knepley 
500ed4e918cSMatthew G Knepley   Level: developer
501ed4e918cSMatthew G Knepley 
502ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
503ed4e918cSMatthew G Knepley @*/
504ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
505341c9bdaSMatthew G Knepley {
506341c9bdaSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
507341c9bdaSMatthew G Knepley   const PetscInt dim = da->dim;
5084e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys/*, mz = da->Ze - da->Zs*/;
509ed4e918cSMatthew 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;
510341c9bdaSMatthew G Knepley 
511341c9bdaSMatthew G Knepley   PetscFunctionBegin;
512ed4e918cSMatthew G Knepley   *cell = -1;
513ed4e918cSMatthew 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);
514ed4e918cSMatthew 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);
515ed4e918cSMatthew 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);
516ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
517ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
518341c9bdaSMatthew G Knepley }
519341c9bdaSMatthew G Knepley 
52057459e9aSMatthew G Knepley #undef __FUNCT__
52157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
52215d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
52357459e9aSMatthew G Knepley {
52457459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
52557459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
52657459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
52757459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
52857459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
52957459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
53057459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
53157459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
53257459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
53357459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
53457459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
53557459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
53657459e9aSMatthew G Knepley   PetscReal         invDet;
53757459e9aSMatthew G Knepley   PetscErrorCode    ierr;
53857459e9aSMatthew G Knepley 
53957459e9aSMatthew G Knepley   PetscFunctionBegin;
54015d2e57cSJed Brown #if defined(PETSC_USE_DEBUG)
54115d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
5424c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
54315d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
54415d2e57cSJed Brown #endif
54515d2e57cSJed Brown   J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
54615d2e57cSJed Brown   J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
54757459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
54857459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
54957459e9aSMatthew G Knepley   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
55057459e9aSMatthew G Knepley   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
55157459e9aSMatthew G Knepley   ierr = PetscLogFlops(30);CHKERRQ(ierr);
55257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
55357459e9aSMatthew G Knepley }
55457459e9aSMatthew G Knepley 
55557459e9aSMatthew G Knepley #undef __FUNCT__
55657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
55757459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
55857459e9aSMatthew G Knepley {
55957459e9aSMatthew G Knepley   DM                 cdm;
56057459e9aSMatthew G Knepley   Vec                coordinates;
56157459e9aSMatthew G Knepley   const PetscScalar *vertices;
56257459e9aSMatthew G Knepley   PetscInt           dim, d, q;
56357459e9aSMatthew G Knepley   PetscErrorCode     ierr;
56457459e9aSMatthew G Knepley 
56557459e9aSMatthew G Knepley   PetscFunctionBegin;
56657459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
56757459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
5686636e97aSMatthew G Knepley   ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
5696636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
57057459e9aSMatthew G Knepley   ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr);
57157459e9aSMatthew G Knepley   for (d = 0; d < dim; ++d) {
57215d2e57cSJed Brown     v0[d] = PetscRealPart(vertices[d]);
57357459e9aSMatthew G Knepley   }
57457459e9aSMatthew G Knepley   switch(dim) {
57557459e9aSMatthew G Knepley   case 2:
57657459e9aSMatthew G Knepley     for (q = 0; q < quad->numQuadPoints; ++q) {
57757459e9aSMatthew G Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
57857459e9aSMatthew G Knepley     }
57957459e9aSMatthew G Knepley     break;
58057459e9aSMatthew G Knepley   default:
58157459e9aSMatthew G Knepley     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
58257459e9aSMatthew G Knepley   }
58357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
58457459e9aSMatthew G Knepley }
585