xref: /petsc/src/dm/impls/da/dageometry.c (revision da80777bb9afa786145217fc25b54e33c3ebf541)
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 */
115201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
116201c01f8SBarry 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;
133*da80777bSKarl Rupp       PetscInt points[9];
134*da80777bSKarl Rupp 
135*da80777bSKarl Rupp       /* Note: initializer list is not computable at compile time, hence we fill the array manually */
136*da80777bSKarl 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;
137aa1993deSMatthew G Knepley 
138aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
139201c01f8SBarry Smith     } else SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
140aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
141aa1993deSMatthew G Knepley     /* Vertex */
142aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
143aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
144aa1993deSMatthew G Knepley     /* X Face */
145201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
146201c01f8SBarry Smith     else if (dim == 2) {
147aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
148aa1993deSMatthew G Knepley       PetscInt f         = p - xfStart;
149*da80777bSKarl Rupp       PetscInt points[3];
150aa1993deSMatthew G Knepley 
151*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
152aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
153aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
154aa1993deSMatthew G Knepley     } else if (dim == 3) {
155aa1993deSMatthew G Knepley       /* 4 vertices */
156aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
157aa1993deSMatthew G Knepley     }
158aa1993deSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
159aa1993deSMatthew G Knepley     /* Y Face */
160201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
161201c01f8SBarry Smith     else if (dim == 2) {
162aa1993deSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
163aa1993deSMatthew G Knepley       PetscInt f         = p - yfStart;
164*da80777bSKarl Rupp       PetscInt points[3];
165aa1993deSMatthew G Knepley 
166*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2]= f+1;
167aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
168aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
169aa1993deSMatthew G Knepley     } else if (dim == 3) {
170aa1993deSMatthew G Knepley       /* 4 vertices */
171aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
172aa1993deSMatthew G Knepley     }
173aa1993deSMatthew G Knepley   } else {
174aa1993deSMatthew G Knepley     /* Z Face */
175201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
176201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
177201c01f8SBarry Smith     else if (dim == 3) {
178aa1993deSMatthew G Knepley       /* 4 vertices */
179aa1993deSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
180aa1993deSMatthew G Knepley     }
181aa1993deSMatthew G Knepley   }
182aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
183aa1993deSMatthew G Knepley }
184aa1993deSMatthew G Knepley 
185aa1993deSMatthew G Knepley #undef __FUNCT__
186aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
187aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
188aa1993deSMatthew G Knepley {
189aa1993deSMatthew G Knepley   PetscErrorCode ierr;
190aa1993deSMatthew G Knepley 
191aa1993deSMatthew G Knepley   PetscFunctionBegin;
192aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
193aa1993deSMatthew G Knepley   PetscValidIntPointer(n,4);
194aa1993deSMatthew G Knepley   PetscValidPointer(closure, 5);
195aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
196aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
197aa1993deSMatthew G Knepley }
198aa1993deSMatthew G Knepley 
199aa1993deSMatthew G Knepley #undef __FUNCT__
200aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
201aa1993deSMatthew G Knepley /* If you did not pass PETSC_NULL for 'values', you must call DMDARestoreClosureScalar() */
202aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
203aa1993deSMatthew G Knepley {
204aa1993deSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
205aa1993deSMatthew G Knepley   PetscInt       dim = da->dim;
20657459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
20795babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
20857459e9aSMatthew G Knepley   PetscErrorCode ierr;
20957459e9aSMatthew G Knepley 
21057459e9aSMatthew G Knepley   PetscFunctionBegin;
21157459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
212aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
21357459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
21457459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
21557459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
21657459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
21757459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
21857459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
21957459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
22057459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
22157459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
22257459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
22357459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
22495babc02SJed Brown   zfStart = yfEnd;
22557459e9aSMatthew 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);
22657459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
22757459e9aSMatthew G Knepley     /* Cell */
228201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
229201c01f8SBarry Smith     else if (dim == 2) {
23057459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
23157459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
23257459e9aSMatthew G Knepley          Bottom y-face same order as cells
23357459e9aSMatthew G Knepley          Left x-face follows same order as cells
23457459e9aSMatthew G Knepley          We number the quad:
23557459e9aSMatthew G Knepley 
23657459e9aSMatthew G Knepley            8--3--7
23757459e9aSMatthew G Knepley            |     |
23857459e9aSMatthew G Knepley            4  0  2
23957459e9aSMatthew G Knepley            |     |
24057459e9aSMatthew G Knepley            5--1--6
24157459e9aSMatthew G Knepley       */
24257459e9aSMatthew G Knepley       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
24357459e9aSMatthew G Knepley       PetscInt v         = cy*nVx + cx +  vStart;
24457459e9aSMatthew G Knepley       PetscInt xf        = cy*nxF + cx + xfStart;
24557459e9aSMatthew G Knepley       PetscInt yf        = c + yfStart;
246*da80777bSKarl Rupp       PetscInt points[9];
24757459e9aSMatthew G Knepley 
248*da80777bSKarl 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;
24957459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr);
25057459e9aSMatthew G Knepley     } else {
25157459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
25257459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
25357459e9aSMatthew G Knepley          Back z-face follows same order as cells
25457459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
25557459e9aSMatthew G Knepley          Left x-face follows same order as cells
25657459e9aSMatthew G Knepley 
25757459e9aSMatthew G Knepley               14-----13
25857459e9aSMatthew G Knepley               /|    /|
25957459e9aSMatthew G Knepley              / | 2 / |
26057459e9aSMatthew G Knepley             / 5|  /  |
26157459e9aSMatthew G Knepley           10-----9  4|
26257459e9aSMatthew G Knepley            |  11-|---12
26357459e9aSMatthew G Knepley            |6 /  |  /
26457459e9aSMatthew G Knepley            | /1 3| /
26557459e9aSMatthew G Knepley            |/    |/
26657459e9aSMatthew G Knepley            7-----8
26757459e9aSMatthew G Knepley       */
26857459e9aSMatthew G Knepley       PetscInt c          = p - cStart;
269*da80777bSKarl Rupp       PetscInt points[15];
270*da80777bSKarl Rupp 
271*da80777bSKarl 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;
272*da80777bSKarl 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;
273*da80777bSKarl Rupp       points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
27457459e9aSMatthew G Knepley 
27557459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
27657459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr);
27757459e9aSMatthew G Knepley     }
27857459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
27957459e9aSMatthew G Knepley     /* Vertex */
28057459e9aSMatthew G Knepley     ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr);
28157459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
28257459e9aSMatthew G Knepley     /* X Face */
283201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
284201c01f8SBarry Smith     else if (dim == 2) {
28557459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
28657459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
287*da80777bSKarl Rupp       PetscInt points[3];
28857459e9aSMatthew G Knepley 
289*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
29057459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
29157459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
29257459e9aSMatthew G Knepley     } else if (dim == 3) {
29357459e9aSMatthew G Knepley       /* 4 vertices */
29457459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
29557459e9aSMatthew G Knepley     }
29657459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
29757459e9aSMatthew G Knepley     /* Y Face */
298201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
299201c01f8SBarry Smith     else if (dim == 2) {
30057459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
30157459e9aSMatthew G Knepley       PetscInt f         = p - yfStart;
302*da80777bSKarl Rupp       PetscInt points[3];
30357459e9aSMatthew G Knepley 
304*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
30557459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
30657459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
30757459e9aSMatthew G Knepley     } else if (dim == 3) {
30857459e9aSMatthew G Knepley       /* 4 vertices */
30957459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
31057459e9aSMatthew G Knepley     }
31157459e9aSMatthew G Knepley   } else {
31257459e9aSMatthew G Knepley     /* Z Face */
313201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
314201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
315201c01f8SBarry Smith     else if (dim == 3) {
31657459e9aSMatthew G Knepley       /* 4 vertices */
31757459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
31857459e9aSMatthew G Knepley     }
31957459e9aSMatthew G Knepley   }
320aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
321aa1993deSMatthew G Knepley }
322aa1993deSMatthew G Knepley 
323aa1993deSMatthew G Knepley #undef __FUNCT__
324aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
325aa1993deSMatthew G Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
326aa1993deSMatthew G Knepley {
327aa1993deSMatthew G Knepley   PetscScalar *vArray;
328aa1993deSMatthew G Knepley   PetscErrorCode ierr;
329aa1993deSMatthew G Knepley 
330aa1993deSMatthew G Knepley   PetscFunctionBegin;
331aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
332aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
333aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
334aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
335aa1993deSMatthew G Knepley   ierr = DMDAGetClosureScalar(dm,section,p,vArray,values);CHKERRQ(ierr);
33657459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
33757459e9aSMatthew G Knepley   PetscFunctionReturn(0);
33857459e9aSMatthew G Knepley }
33957459e9aSMatthew G Knepley 
34057459e9aSMatthew G Knepley #undef __FUNCT__
341aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
342aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
343aa1993deSMatthew G Knepley {
344aa1993deSMatthew G Knepley   PetscErrorCode ierr;
345aa1993deSMatthew G Knepley   PetscInt count;
346aa1993deSMatthew G Knepley 
347aa1993deSMatthew G Knepley   PetscFunctionBegin;
348aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
349aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
350aa1993deSMatthew G Knepley   count = 0;                    /* We are lying about the count */
351854432f1SMatthew G Knepley   ierr = DMRestoreWorkArray(dm, count, PETSC_SCALAR, (void *) values);CHKERRQ(ierr);
352aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
353aa1993deSMatthew G Knepley }
354aa1993deSMatthew G Knepley 
355aa1993deSMatthew G Knepley #undef __FUNCT__
356aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
357aa1993deSMatthew G Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
358aa1993deSMatthew G Knepley {
359aa1993deSMatthew G Knepley   PetscErrorCode ierr;
360aa1993deSMatthew G Knepley 
361aa1993deSMatthew G Knepley   PetscFunctionBegin;
362aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
363aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
364aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
365aa1993deSMatthew G Knepley   ierr = DMDARestoreClosureScalar(dm,section,p,PETSC_NULL,values);CHKERRQ(ierr);
366aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
367aa1993deSMatthew G Knepley }
368aa1993deSMatthew G Knepley 
369aa1993deSMatthew G Knepley #undef __FUNCT__
370aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
371aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
37257459e9aSMatthew G Knepley {
37357459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
37457459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
37557459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
37695babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
37757459e9aSMatthew G Knepley   PetscErrorCode ierr;
37857459e9aSMatthew G Knepley 
37957459e9aSMatthew G Knepley   PetscFunctionBegin;
38057459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
381aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
38257459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
38357459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
38457459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
38557459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
38657459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
38757459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
38857459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
38957459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
39057459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
39157459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
39257459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
39395babc02SJed Brown   zfStart = yfEnd;
39457459e9aSMatthew 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);
39557459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
39657459e9aSMatthew G Knepley     /* Cell */
397201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
398201c01f8SBarry Smith     else if (dim == 2) {
39957459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
40057459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
40157459e9aSMatthew G Knepley          Bottom y-face same order as cells
40257459e9aSMatthew G Knepley          Left x-face follows same order as cells
40357459e9aSMatthew G Knepley          We number the quad:
40457459e9aSMatthew G Knepley 
40557459e9aSMatthew G Knepley            8--3--7
40657459e9aSMatthew G Knepley            |     |
40757459e9aSMatthew G Knepley            4  0  2
40857459e9aSMatthew G Knepley            |     |
40957459e9aSMatthew G Knepley            5--1--6
41057459e9aSMatthew G Knepley       */
41157459e9aSMatthew G Knepley       PetscInt c         = p - cStart;
412*da80777bSKarl Rupp       PetscInt points[9];
41357459e9aSMatthew G Knepley 
414*da80777bSKarl 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;
415*da80777bSKarl Rupp       points[6] = c+vStart+1; points[7] = c+vStart+nVx+1; points[8] = c+vStart+nVx+0;
41657459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
41757459e9aSMatthew G Knepley     } else {
41857459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
41957459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
42057459e9aSMatthew G Knepley          Back z-face follows same order as cells
42157459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
42257459e9aSMatthew G Knepley          Left x-face follows same order as cells
42357459e9aSMatthew G Knepley 
42457459e9aSMatthew G Knepley               14-----13
42557459e9aSMatthew G Knepley               /|    /|
42657459e9aSMatthew G Knepley              / | 2 / |
42757459e9aSMatthew G Knepley             / 5|  /  |
42857459e9aSMatthew G Knepley           10-----9  4|
42957459e9aSMatthew G Knepley            |  11-|---12
43057459e9aSMatthew G Knepley            |6 /  |  /
43157459e9aSMatthew G Knepley            | /1 3| /
43257459e9aSMatthew G Knepley            |/    |/
43357459e9aSMatthew G Knepley            7-----8
43457459e9aSMatthew G Knepley       */
43557459e9aSMatthew G Knepley       PetscInt c          = p - cStart;
436*da80777bSKarl Rupp       PetscInt points[15];
43757459e9aSMatthew G Knepley 
438*da80777bSKarl 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;
439*da80777bSKarl 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;
440*da80777bSKarl Rupp       points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
44157459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
44257459e9aSMatthew G Knepley     }
44357459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
44457459e9aSMatthew G Knepley     /* Vertex */
44557459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
44657459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
44757459e9aSMatthew G Knepley     /* X Face */
448201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
449201c01f8SBarry Smith     else if (dim == 2) {
45057459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
45157459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
452*da80777bSKarl Rupp       PetscInt points[3];
45357459e9aSMatthew G Knepley 
454*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
45557459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
45657459e9aSMatthew G Knepley     } else if (dim == 3) {
45757459e9aSMatthew G Knepley       /* 4 vertices */
45857459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
45957459e9aSMatthew G Knepley     }
46057459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
46157459e9aSMatthew G Knepley     /* Y Face */
462201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
463201c01f8SBarry Smith     else if (dim == 2) {
46457459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
46557459e9aSMatthew G Knepley       PetscInt f         = p - yfStart;
466*da80777bSKarl Rupp       PetscInt points[3];
46757459e9aSMatthew G Knepley 
468*da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
46957459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
47057459e9aSMatthew G Knepley     } else if (dim == 3) {
47157459e9aSMatthew G Knepley       /* 4 vertices */
47257459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
47357459e9aSMatthew G Knepley     }
47457459e9aSMatthew G Knepley   } else {
47557459e9aSMatthew G Knepley     /* Z Face */
476201c01f8SBarry Smith     if (dim == 1) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
477201c01f8SBarry Smith     else if (dim == 2) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
478201c01f8SBarry Smith     else if (dim == 3) {
47957459e9aSMatthew G Knepley       /* 4 vertices */
48057459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
48157459e9aSMatthew G Knepley     }
48257459e9aSMatthew G Knepley   }
483aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
484aa1993deSMatthew G Knepley }
485aa1993deSMatthew G Knepley 
486aa1993deSMatthew G Knepley #undef __FUNCT__
487aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
488aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
489aa1993deSMatthew G Knepley {
490aa1993deSMatthew G Knepley   PetscScalar   *vArray;
491aa1993deSMatthew G Knepley   PetscErrorCode ierr;
492aa1993deSMatthew G Knepley 
493aa1993deSMatthew G Knepley   PetscFunctionBegin;
494aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
495aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
496aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
497aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
498aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
49957459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
50057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
50157459e9aSMatthew G Knepley }
50257459e9aSMatthew G Knepley 
503ed4e918cSMatthew G Knepley #undef __FUNCT__
504ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
505ed4e918cSMatthew G Knepley /*@
506f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
507341c9bdaSMatthew G Knepley 
508ed4e918cSMatthew G Knepley   Not Collective
509ed4e918cSMatthew G Knepley 
510ed4e918cSMatthew G Knepley   Input Parameter:
511ed4e918cSMatthew G Knepley + da - the distributed array
512ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
513ed4e918cSMatthew G Knepley 
514ed4e918cSMatthew G Knepley   Output Parameter:
515ed4e918cSMatthew G Knepley . cell - the local cell number
516ed4e918cSMatthew G Knepley 
517ed4e918cSMatthew G Knepley   Level: developer
518ed4e918cSMatthew G Knepley 
519ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
520ed4e918cSMatthew G Knepley @*/
521ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
522341c9bdaSMatthew G Knepley {
523341c9bdaSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
524341c9bdaSMatthew G Knepley   const PetscInt dim = da->dim;
5254e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys/*, mz = da->Ze - da->Zs*/;
526ed4e918cSMatthew 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;
527341c9bdaSMatthew G Knepley 
528341c9bdaSMatthew G Knepley   PetscFunctionBegin;
529ed4e918cSMatthew G Knepley   *cell = -1;
530ed4e918cSMatthew 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);
531ed4e918cSMatthew 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);
532ed4e918cSMatthew 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);
533ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
534ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
535341c9bdaSMatthew G Knepley }
536341c9bdaSMatthew G Knepley 
53757459e9aSMatthew G Knepley #undef __FUNCT__
53857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
53915d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
54057459e9aSMatthew G Knepley {
54157459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
54257459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
54357459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
54457459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
54557459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
54657459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
54757459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
54857459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
54957459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
55057459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
55157459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
55257459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
55357459e9aSMatthew G Knepley   PetscReal         invDet;
55457459e9aSMatthew G Knepley   PetscErrorCode    ierr;
55557459e9aSMatthew G Knepley 
55657459e9aSMatthew G Knepley   PetscFunctionBegin;
55715d2e57cSJed Brown #if defined(PETSC_USE_DEBUG)
55815d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
5594c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
56015d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
56115d2e57cSJed Brown #endif
56215d2e57cSJed Brown   J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
56315d2e57cSJed Brown   J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
56457459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
56557459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
56657459e9aSMatthew G Knepley   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
56757459e9aSMatthew G Knepley   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
56857459e9aSMatthew G Knepley   ierr = PetscLogFlops(30);CHKERRQ(ierr);
56957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
57057459e9aSMatthew G Knepley }
57157459e9aSMatthew G Knepley 
57257459e9aSMatthew G Knepley #undef __FUNCT__
57357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
57457459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
57557459e9aSMatthew G Knepley {
57657459e9aSMatthew G Knepley   DM                 cdm;
57757459e9aSMatthew G Knepley   Vec                coordinates;
57857459e9aSMatthew G Knepley   const PetscScalar *vertices;
57957459e9aSMatthew G Knepley   PetscInt           dim, d, q;
58057459e9aSMatthew G Knepley   PetscErrorCode     ierr;
58157459e9aSMatthew G Knepley 
58257459e9aSMatthew G Knepley   PetscFunctionBegin;
58357459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
58457459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
5856636e97aSMatthew G Knepley   ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
5866636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
58757459e9aSMatthew G Knepley   ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr);
58857459e9aSMatthew G Knepley   for (d = 0; d < dim; ++d) {
58915d2e57cSJed Brown     v0[d] = PetscRealPart(vertices[d]);
59057459e9aSMatthew G Knepley   }
59157459e9aSMatthew G Knepley   switch (dim) {
59257459e9aSMatthew G Knepley   case 2:
59357459e9aSMatthew G Knepley     for (q = 0; q < quad->numQuadPoints; ++q) {
59457459e9aSMatthew G Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
59557459e9aSMatthew G Knepley     }
59657459e9aSMatthew G Knepley     break;
59757459e9aSMatthew G Knepley   default:
59857459e9aSMatthew G Knepley     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
59957459e9aSMatthew G Knepley   }
60057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
60157459e9aSMatthew G Knepley }
602