xref: /petsc/src/dm/impls/da/dageometry.c (revision c73cfb54f34510b9dcefc99e398efa3b88d5dde5)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>     /*I  "petscdmda.h"   I*/
257459e9aSMatthew G Knepley 
357459e9aSMatthew G Knepley #undef __FUNCT__
40a3ada39SMatthew G. Knepley #define __FUNCT__ "FillClosureArray_Static"
5d7a12f93SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Static(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, PetscInt *csize, PetscScalar **array)
657459e9aSMatthew G Knepley {
757459e9aSMatthew G Knepley   PetscScalar    *a;
80a3ada39SMatthew G. Knepley   PetscInt       pStart, pEnd, size = 0, dof, off, d, k, i;
957459e9aSMatthew G Knepley   PetscErrorCode ierr;
1057459e9aSMatthew G Knepley 
1157459e9aSMatthew G Knepley   PetscFunctionBegin;
120a3ada39SMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
1357459e9aSMatthew G Knepley   for (i = 0; i < nP; ++i) {
140a3ada39SMatthew G. Knepley     const PetscInt p = points[i];
150a3ada39SMatthew G. Knepley 
160a3ada39SMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
170a3ada39SMatthew G. Knepley     ierr  = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1857459e9aSMatthew G Knepley     size += dof;
1957459e9aSMatthew G Knepley   }
200a3ada39SMatthew G. Knepley   if (csize) *csize = size;
210a3ada39SMatthew G. Knepley   if (array) {
22aa1993deSMatthew G Knepley     ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr);
2357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
240a3ada39SMatthew G. Knepley       const PetscInt p = points[i];
250a3ada39SMatthew G. Knepley 
260a3ada39SMatthew G. Knepley       if ((p < pStart) || (p >= pEnd)) continue;
270a3ada39SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
280a3ada39SMatthew G. Knepley       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
2957459e9aSMatthew G Knepley 
308865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d];
3157459e9aSMatthew G Knepley     }
3257459e9aSMatthew G Knepley     *array = a;
330a3ada39SMatthew G. Knepley   }
3457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
3557459e9aSMatthew G Knepley }
3657459e9aSMatthew G Knepley 
3757459e9aSMatthew G Knepley #undef __FUNCT__
3857459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
3957459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
4057459e9aSMatthew G Knepley {
4157459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
4257459e9aSMatthew G Knepley   PetscErrorCode ierr;
4357459e9aSMatthew G Knepley 
4457459e9aSMatthew G Knepley   PetscFunctionBegin;
4557459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
4657459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
4757459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4857459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
4957459e9aSMatthew G Knepley 
508865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k];
5157459e9aSMatthew G Knepley     }
5257459e9aSMatthew G Knepley   } else {
5357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
5457459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
5557459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
5657459e9aSMatthew G Knepley 
578865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k];
5857459e9aSMatthew G Knepley     }
5957459e9aSMatthew G Knepley   }
6057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
6157459e9aSMatthew G Knepley }
6257459e9aSMatthew G Knepley 
6357459e9aSMatthew G Knepley #undef __FUNCT__
64aa1993deSMatthew G Knepley #define __FUNCT__ "GetPointArray_Private"
65aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
66aa1993deSMatthew G Knepley {
67aa1993deSMatthew G Knepley   PetscErrorCode ierr;
68aa1993deSMatthew G Knepley   PetscInt       *work;
69aa1993deSMatthew G Knepley 
70aa1993deSMatthew G Knepley   PetscFunctionBegin;
71aa1993deSMatthew G Knepley   if (rn) *rn = n;
72aa1993deSMatthew G Knepley   if (rpoints) {
73aa1993deSMatthew G Knepley     ierr     = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
74aa1993deSMatthew G Knepley     ierr     = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
75aa1993deSMatthew G Knepley     *rpoints = work;
76aa1993deSMatthew G Knepley   }
77aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
78aa1993deSMatthew G Knepley }
79aa1993deSMatthew G Knepley 
80aa1993deSMatthew G Knepley #undef __FUNCT__
81aa1993deSMatthew G Knepley #define __FUNCT__ "RestorePointArray_Private"
82aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
83aa1993deSMatthew G Knepley {
84aa1993deSMatthew G Knepley   PetscErrorCode ierr;
85aa1993deSMatthew G Knepley 
86aa1993deSMatthew G Knepley   PetscFunctionBegin;
87aa1993deSMatthew G Knepley   if (rn) *rn = 0;
888865f1eaSKarl Rupp   if (rpoints) {
898865f1eaSKarl Rupp     ierr = DMRestoreWorkArray(dm,*rn, PETSC_INT, (void*) rpoints);CHKERRQ(ierr);
908865f1eaSKarl Rupp   }
91aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
92aa1993deSMatthew G Knepley }
93aa1993deSMatthew G Knepley 
94aa1993deSMatthew G Knepley #undef __FUNCT__
956693a731SMatthew G. Knepley #define __FUNCT__ "DMDAGetTransitiveClosure"
966693a731SMatthew G. Knepley PetscErrorCode DMDAGetTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure)
976693a731SMatthew G. Knepley {
98*c73cfb54SMatthew G. Knepley   PetscInt       dim = dm->dim;
996693a731SMatthew G. Knepley   PetscInt       nVx, nVy, nVz, nxF, nXF, nyF, nYF, nzF, nZF;
1006693a731SMatthew G. Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
1016693a731SMatthew G. Knepley   PetscErrorCode ierr;
1026693a731SMatthew G. Knepley 
1036693a731SMatthew G. Knepley   PetscFunctionBegin;
1046693a731SMatthew G. Knepley   if (!useClosure) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Star operation is not yet supported");
1056693a731SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1066693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
1076693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm,  0,  &cStart, &cEnd);CHKERRQ(ierr);
1086693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm,  1,  &fStart, &fEnd);CHKERRQ(ierr);
1096693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
1106693a731SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, NULL);CHKERRQ(ierr);
1116693a731SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
1126693a731SMatthew G. Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
1136693a731SMatthew G. Knepley   yfStart = xfEnd;
1146693a731SMatthew G. Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
1156693a731SMatthew G. Knepley   if ((p >= cStart) || (p < cEnd)) {
1166693a731SMatthew G. Knepley     /* Cell */
1176693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
1186693a731SMatthew G. Knepley     else if (dim == 2) {
1196693a731SMatthew G. Knepley       /* 4 faces, 4 vertices
1206693a731SMatthew G. Knepley          Bottom-left vertex follows same order as cells
1216693a731SMatthew G. Knepley          Bottom y-face same order as cells
1226693a731SMatthew G. Knepley          Left x-face follows same order as cells
1236693a731SMatthew G. Knepley          We number the quad:
1246693a731SMatthew G. Knepley 
1256693a731SMatthew G. Knepley            8--3--7
1266693a731SMatthew G. Knepley            |     |
1276693a731SMatthew G. Knepley            4  0  2
1286693a731SMatthew G. Knepley            |     |
1296693a731SMatthew G. Knepley            5--1--6
1306693a731SMatthew G. Knepley       */
1316693a731SMatthew G. Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
1326693a731SMatthew G. Knepley       PetscInt v  = cy*nVx + cx +  vStart;
1336693a731SMatthew G. Knepley       PetscInt xf = cy*nxF + cx + xfStart;
1346693a731SMatthew G. Knepley       PetscInt yf = c + yfStart;
1356693a731SMatthew G. Knepley 
1366693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 9;}
1376693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
1386693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = yf; (*closure)[2] = xf+1; (*closure)[3] = yf+nyF; (*closure)[4] = xf+0; (*closure)[5] = v+0; (*closure)[6]= v+1; (*closure)[7] = v+nVx+1; (*closure)[8] = v+nVx+0;
1396693a731SMatthew G. Knepley     } else {
1406693a731SMatthew G. Knepley       /* 6 faces, 12 edges, 8 vertices
1416693a731SMatthew G. Knepley          Bottom-left vertex follows same order as cells
1426693a731SMatthew G. Knepley          Bottom y-face same order as cells
1436693a731SMatthew G. Knepley          Left x-face follows same order as cells
1446693a731SMatthew G. Knepley          We number the quad:
1456693a731SMatthew G. Knepley 
1466693a731SMatthew G. Knepley            8--3--7
1476693a731SMatthew G. Knepley            |     |
1486693a731SMatthew G. Knepley            4  0  2
1496693a731SMatthew G. Knepley            |     |
1506693a731SMatthew G. Knepley            5--1--6
1516693a731SMatthew G. Knepley       */
152be8e0784SMatthew G. Knepley #if 0
1536693a731SMatthew G. Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1), cz = c / ((nVx-1)*(nVy-1));
1546693a731SMatthew G. Knepley       PetscInt v  = cy*nVx + cx +  vStart;
1556693a731SMatthew G. Knepley       PetscInt xf = cy*nxF + cx + xfStart;
1566693a731SMatthew G. Knepley       PetscInt yf = c + yfStart;
1576693a731SMatthew G. Knepley 
1586693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 26;}
1596693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
160be8e0784SMatthew G. Knepley #endif
161be8e0784SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
1626693a731SMatthew G. Knepley     }
1636693a731SMatthew G. Knepley   } else if ((p >= vStart) || (p < vEnd)) {
1646693a731SMatthew G. Knepley     /* Vertex */
1656693a731SMatthew G. Knepley     if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 1;}
1666693a731SMatthew G. Knepley     if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
1676693a731SMatthew G. Knepley     (*closure)[0] = p;
1686693a731SMatthew G. Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
1696693a731SMatthew G. Knepley     /* X Face */
1706693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
1716693a731SMatthew G. Knepley     else if (dim == 2) {
1726693a731SMatthew G. Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
1736693a731SMatthew G. Knepley       PetscInt f = p - xfStart;
1746693a731SMatthew G. Knepley 
1756693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;}
1766693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
1776693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = f; (*closure)[2] = f+nVx;
1786693a731SMatthew G. Knepley     } else if (dim == 3) {
1796693a731SMatthew G. Knepley       /* 4 vertices */
1806693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
1816693a731SMatthew G. Knepley     }
1826693a731SMatthew G. Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
1836693a731SMatthew G. Knepley     /* Y Face */
1846693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
1856693a731SMatthew G. Knepley     else if (dim == 2) {
1866693a731SMatthew G. Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
1876693a731SMatthew G. Knepley       PetscInt f = p - yfStart;
1886693a731SMatthew G. Knepley 
1896693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;}
1906693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
1916693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = f; (*closure)[2]= f+1;
1926693a731SMatthew G. Knepley     } else if (dim == 3) {
1936693a731SMatthew G. Knepley       /* 4 vertices */
1946693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
1956693a731SMatthew G. Knepley     }
1966693a731SMatthew G. Knepley   } else {
1976693a731SMatthew G. Knepley     /* Z Face */
1986693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
1996693a731SMatthew G. Knepley     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
2006693a731SMatthew G. Knepley     else if (dim == 3) {
2016693a731SMatthew G. Knepley       /* 4 vertices */
2026693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
2036693a731SMatthew G. Knepley     }
2046693a731SMatthew G. Knepley   }
2056693a731SMatthew G. Knepley   PetscFunctionReturn(0);
2066693a731SMatthew G. Knepley }
2076693a731SMatthew G. Knepley 
2086693a731SMatthew G. Knepley #undef __FUNCT__
2096693a731SMatthew G. Knepley #define __FUNCT__ "DMDARestoreTransitiveClosure"
2106693a731SMatthew G. Knepley PetscErrorCode DMDARestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure)
2116693a731SMatthew G. Knepley {
2126693a731SMatthew G. Knepley   PetscErrorCode ierr;
2136693a731SMatthew G. Knepley 
2146693a731SMatthew G. Knepley   PetscFunctionBegin;
2156693a731SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, closure);CHKERRQ(ierr);
2166693a731SMatthew G. Knepley   PetscFunctionReturn(0);
2176693a731SMatthew G. Knepley }
2186693a731SMatthew G. Knepley 
2196693a731SMatthew G. Knepley #undef __FUNCT__
220aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosure"
221aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
22257459e9aSMatthew G Knepley {
223*c73cfb54SMatthew G. Knepley   PetscInt       dim = dm->dim;
224aa1993deSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
225aa1993deSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
226aa1993deSMatthew G Knepley   PetscErrorCode ierr;
227aa1993deSMatthew G Knepley 
228aa1993deSMatthew G Knepley   PetscFunctionBegin;
229aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
230aa1993deSMatthew G Knepley   if (n) PetscValidIntPointer(n,4);
231aa1993deSMatthew G Knepley   if (closure) PetscValidPointer(closure, 5);
232aa1993deSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
23382f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
234aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
235aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
236aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
237aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
2380298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
239aa1993deSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
240aa1993deSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
241aa1993deSMatthew G Knepley   yfStart = xfEnd;
24282f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
243aa1993deSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
244aa1993deSMatthew G Knepley     /* Cell */
24582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
246201c01f8SBarry Smith     else if (dim == 2) {
247aa1993deSMatthew G Knepley       /* 4 faces, 4 vertices
248aa1993deSMatthew G Knepley          Bottom-left vertex follows same order as cells
249aa1993deSMatthew G Knepley          Bottom y-face same order as cells
250aa1993deSMatthew G Knepley          Left x-face follows same order as cells
251aa1993deSMatthew G Knepley          We number the quad:
252aa1993deSMatthew G Knepley 
253aa1993deSMatthew G Knepley            8--3--7
254aa1993deSMatthew G Knepley            |     |
255aa1993deSMatthew G Knepley            4  0  2
256aa1993deSMatthew G Knepley            |     |
257aa1993deSMatthew G Knepley            5--1--6
258aa1993deSMatthew G Knepley       */
259aa1993deSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
260aa1993deSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
261aa1993deSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
262aa1993deSMatthew G Knepley       PetscInt yf = c + yfStart;
263da80777bSKarl Rupp       PetscInt points[9];
264da80777bSKarl Rupp 
265da80777bSKarl Rupp       /* Note: initializer list is not computable at compile time, hence we fill the array manually */
266da80777bSKarl 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;
267aa1993deSMatthew G Knepley 
268aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
26982f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
270aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
271aa1993deSMatthew G Knepley     /* Vertex */
272aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
273aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
274aa1993deSMatthew G Knepley     /* X Face */
27582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
276201c01f8SBarry Smith     else if (dim == 2) {
277aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
278aa1993deSMatthew G Knepley       PetscInt f = p - xfStart;
279da80777bSKarl Rupp       PetscInt points[3];
280aa1993deSMatthew G Knepley 
281da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
28282f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
283aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
284aa1993deSMatthew G Knepley     } else if (dim == 3) {
285aa1993deSMatthew G Knepley       /* 4 vertices */
28682f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
287aa1993deSMatthew G Knepley     }
288aa1993deSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
289aa1993deSMatthew G Knepley     /* Y Face */
29082f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
291201c01f8SBarry Smith     else if (dim == 2) {
292aa1993deSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
293aa1993deSMatthew G Knepley       PetscInt f = p - yfStart;
294da80777bSKarl Rupp       PetscInt points[3];
295aa1993deSMatthew G Knepley 
296da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2]= f+1;
29782f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
298aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
299aa1993deSMatthew G Knepley     } else if (dim == 3) {
300aa1993deSMatthew G Knepley       /* 4 vertices */
30182f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
302aa1993deSMatthew G Knepley     }
303aa1993deSMatthew G Knepley   } else {
304aa1993deSMatthew G Knepley     /* Z Face */
30582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
30682f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
307201c01f8SBarry Smith     else if (dim == 3) {
308aa1993deSMatthew G Knepley       /* 4 vertices */
30982f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
310aa1993deSMatthew G Knepley     }
311aa1993deSMatthew G Knepley   }
312aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
313aa1993deSMatthew G Knepley }
314aa1993deSMatthew G Knepley 
315aa1993deSMatthew G Knepley #undef __FUNCT__
316aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
317e5c84f05SJed Brown /* Zeros n and closure. */
318aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
319aa1993deSMatthew G Knepley {
320aa1993deSMatthew G Knepley   PetscErrorCode ierr;
321aa1993deSMatthew G Knepley 
322aa1993deSMatthew G Knepley   PetscFunctionBegin;
323aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
324e5c84f05SJed Brown   if (n) PetscValidIntPointer(n,4);
325e5c84f05SJed Brown   if (closure) PetscValidPointer(closure, 5);
326aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
327aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
328aa1993deSMatthew G Knepley }
329aa1993deSMatthew G Knepley 
330aa1993deSMatthew G Knepley #undef __FUNCT__
331aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
332d7a12f93SMatthew G. Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, PetscScalar **values)
333aa1993deSMatthew G Knepley {
334*c73cfb54SMatthew G. Knepley   PetscInt       dim = dm->dim;
33557459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
33695babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
33757459e9aSMatthew G Knepley   PetscErrorCode ierr;
33857459e9aSMatthew G Knepley 
33957459e9aSMatthew G Knepley   PetscFunctionBegin;
34057459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
341aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
3420a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
34357459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
34482f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
34557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
34657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
34757459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
34857459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
3490298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
35057459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
35157459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
35257459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
35395babc02SJed Brown   zfStart = yfEnd;
35482f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
35557459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
35657459e9aSMatthew G Knepley     /* Cell */
35782f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
358201c01f8SBarry Smith     else if (dim == 2) {
35957459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
36057459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
36157459e9aSMatthew G Knepley          Bottom y-face same order as cells
36257459e9aSMatthew G Knepley          Left x-face follows same order as cells
36357459e9aSMatthew G Knepley          We number the quad:
36457459e9aSMatthew G Knepley 
36557459e9aSMatthew G Knepley            8--3--7
36657459e9aSMatthew G Knepley            |     |
36757459e9aSMatthew G Knepley            4  0  2
36857459e9aSMatthew G Knepley            |     |
36957459e9aSMatthew G Knepley            5--1--6
37057459e9aSMatthew G Knepley       */
37157459e9aSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
37257459e9aSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
373312b75c6SMatthew G. Knepley       PetscInt xf = cx*nxF + cy + xfStart;
37457459e9aSMatthew G Knepley       PetscInt yf = c + yfStart;
375da80777bSKarl Rupp       PetscInt points[9];
37657459e9aSMatthew G Knepley 
377312b75c6SMatthew G. Knepley       points[0] = p;
378312b75c6SMatthew G. Knepley       points[1] = yf;  points[2] = xf+nxF; points[3] = yf+nyF;  points[4] = xf;
379312b75c6SMatthew G. Knepley       points[5] = v+0; points[6] = v+1;    points[7] = v+nVx+1; points[8] = v+nVx+0;
3800a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 9, points, vArray, csize, values);CHKERRQ(ierr);
38157459e9aSMatthew G Knepley     } else {
38257459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
38357459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
38457459e9aSMatthew G Knepley          Back z-face follows same order as cells
38557459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
38657459e9aSMatthew G Knepley          Left x-face follows same order as cells
38757459e9aSMatthew G Knepley 
38857459e9aSMatthew G Knepley               14-----13
38957459e9aSMatthew G Knepley               /|    /|
39057459e9aSMatthew G Knepley              / | 2 / |
39157459e9aSMatthew G Knepley             / 5|  /  |
39257459e9aSMatthew G Knepley           10-----9  4|
39357459e9aSMatthew G Knepley            |  11-|---12
39457459e9aSMatthew G Knepley            |6 /  |  /
39557459e9aSMatthew G Knepley            | /1 3| /
39657459e9aSMatthew G Knepley            |/    |/
39757459e9aSMatthew G Knepley            7-----8
39857459e9aSMatthew G Knepley       */
39957459e9aSMatthew G Knepley       PetscInt c = p - cStart;
400da80777bSKarl Rupp       PetscInt points[15];
401da80777bSKarl Rupp 
402da80777bSKarl 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;
403da80777bSKarl 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;
404da80777bSKarl Rupp       points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
40557459e9aSMatthew G Knepley 
40682f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4070a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 15, points, vArray, csize, values);CHKERRQ(ierr);
40857459e9aSMatthew G Knepley     }
40957459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
41057459e9aSMatthew G Knepley     /* Vertex */
4110a3ada39SMatthew G. Knepley     ierr = FillClosureArray_Static(dm, section, 1, &p, vArray, csize, values);CHKERRQ(ierr);
41257459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
41357459e9aSMatthew G Knepley     /* X Face */
41482f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
415201c01f8SBarry Smith     else if (dim == 2) {
41657459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
41757459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
418da80777bSKarl Rupp       PetscInt points[3];
41957459e9aSMatthew G Knepley 
420da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
42182f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4220a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
42382f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
42457459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
42557459e9aSMatthew G Knepley     /* Y Face */
42682f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
427201c01f8SBarry Smith     else if (dim == 2) {
42857459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
42957459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
430da80777bSKarl Rupp       PetscInt points[3];
43157459e9aSMatthew G Knepley 
432da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
43382f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4340a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
43582f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
43657459e9aSMatthew G Knepley   } else {
43757459e9aSMatthew G Knepley     /* Z Face */
43882f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
43982f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
44082f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
44157459e9aSMatthew G Knepley   }
442aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
443aa1993deSMatthew G Knepley }
444aa1993deSMatthew G Knepley 
445aa1993deSMatthew G Knepley #undef __FUNCT__
446aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
447d7a12f93SMatthew G. Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, PetscScalar **values)
448aa1993deSMatthew G Knepley {
449aa1993deSMatthew G Knepley   PetscScalar    *vArray;
450aa1993deSMatthew G Knepley   PetscErrorCode ierr;
451aa1993deSMatthew G Knepley 
452aa1993deSMatthew G Knepley   PetscFunctionBegin;
453aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
454aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
4550a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
456aa1993deSMatthew G Knepley   ierr = VecGetArray(v, &vArray);CHKERRQ(ierr);
4570a3ada39SMatthew G. Knepley   ierr = DMDAGetClosureScalar(dm, section, p, vArray, csize, values);CHKERRQ(ierr);
45857459e9aSMatthew G Knepley   ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr);
45957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
46057459e9aSMatthew G Knepley }
46157459e9aSMatthew G Knepley 
46257459e9aSMatthew G Knepley #undef __FUNCT__
463aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
464d7a12f93SMatthew G. Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, PetscScalar **values)
465aa1993deSMatthew G Knepley {
466aa1993deSMatthew G Knepley   PetscErrorCode ierr;
467aa1993deSMatthew G Knepley 
468aa1993deSMatthew G Knepley   PetscFunctionBegin;
469aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4700a3ada39SMatthew G. Knepley   PetscValidPointer(values, 6);
4710a3ada39SMatthew G. Knepley   ierr  = DMRestoreWorkArray(dm, csize ? *csize : 0, PETSC_SCALAR, (void*) values);CHKERRQ(ierr);
472aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
473aa1993deSMatthew G Knepley }
474aa1993deSMatthew G Knepley 
475aa1993deSMatthew G Knepley #undef __FUNCT__
476aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
477d7a12f93SMatthew G. Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, PetscScalar **values)
478aa1993deSMatthew G Knepley {
479aa1993deSMatthew G Knepley   PetscErrorCode ierr;
480aa1993deSMatthew G Knepley 
481aa1993deSMatthew G Knepley   PetscFunctionBegin;
482aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
483aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
484aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
4850a3ada39SMatthew G. Knepley   ierr = DMDARestoreClosureScalar(dm, section, p, NULL, csize, values);CHKERRQ(ierr);
486aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
487aa1993deSMatthew G Knepley }
488aa1993deSMatthew G Knepley 
489aa1993deSMatthew G Knepley #undef __FUNCT__
490aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
491aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
49257459e9aSMatthew G Knepley {
493*c73cfb54SMatthew G. Knepley   PetscInt       dim = dm->dim;
4946a7fb054SMatthew G. Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF, nCx, nCy;
49595babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
49657459e9aSMatthew G Knepley   PetscErrorCode ierr;
49757459e9aSMatthew G Knepley 
49857459e9aSMatthew G Knepley   PetscFunctionBegin;
49957459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
500aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
50157459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
50257459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
50382f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
50457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
50557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
50657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
50757459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
5086a7fb054SMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, &nCx, &nCy, NULL, NULL);CHKERRQ(ierr);
5090298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
51057459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
51157459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
51257459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
51395babc02SJed Brown   zfStart = yfEnd;
51482f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
51557459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
51657459e9aSMatthew G Knepley     /* Cell */
51782f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
518201c01f8SBarry Smith     else if (dim == 2) {
51957459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
52057459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
52157459e9aSMatthew G Knepley          Bottom y-face same order as cells
52257459e9aSMatthew G Knepley          Left x-face follows same order as cells
52357459e9aSMatthew G Knepley          We number the quad:
52457459e9aSMatthew G Knepley 
52557459e9aSMatthew G Knepley            8--3--7
52657459e9aSMatthew G Knepley            |     |
52757459e9aSMatthew G Knepley            4  0  2
52857459e9aSMatthew G Knepley            |     |
52957459e9aSMatthew G Knepley            5--1--6
53057459e9aSMatthew G Knepley       */
531312b75c6SMatthew G. Knepley       PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
532312b75c6SMatthew G. Knepley       PetscInt v  = cy*nVx + cx +  vStart;
533312b75c6SMatthew G. Knepley       PetscInt xf = cx*nxF + cy + xfStart;
534312b75c6SMatthew G. Knepley       PetscInt yf = c + yfStart;
535da80777bSKarl Rupp       PetscInt points[9];
53657459e9aSMatthew G Knepley 
5376a7fb054SMatthew G. Knepley       points[0] = p;
538312b75c6SMatthew G. Knepley       points[1] = yf;  points[2] = xf+nxF; points[3] = yf+nyF;  points[4] = xf;
539312b75c6SMatthew G. Knepley       points[5] = v+0; points[6] = v+1;    points[7] = v+nVx+1; points[8] = v+nVx+0;
54057459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
54157459e9aSMatthew G Knepley     } else {
54257459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
54357459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
54457459e9aSMatthew G Knepley          Back z-face follows same order as cells
54557459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
54657459e9aSMatthew G Knepley          Left x-face follows same order as cells
54757459e9aSMatthew G Knepley 
54857459e9aSMatthew G Knepley               14-----13
54957459e9aSMatthew G Knepley               /|    /|
55057459e9aSMatthew G Knepley              / | 2 / |
55157459e9aSMatthew G Knepley             / 5|  /  |
55257459e9aSMatthew G Knepley           10-----9  4|
55357459e9aSMatthew G Knepley            |  11-|---12
55457459e9aSMatthew G Knepley            |6 /  |  /
55557459e9aSMatthew G Knepley            | /1 3| /
55657459e9aSMatthew G Knepley            |/    |/
55757459e9aSMatthew G Knepley            7-----8
55857459e9aSMatthew G Knepley       */
55957459e9aSMatthew G Knepley       PetscInt c = p - cStart;
560da80777bSKarl Rupp       PetscInt points[15];
56157459e9aSMatthew G Knepley 
562da80777bSKarl 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;
563da80777bSKarl 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;
564da80777bSKarl Rupp       points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
56557459e9aSMatthew G Knepley       ierr       = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
56657459e9aSMatthew G Knepley     }
56757459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
56857459e9aSMatthew G Knepley     /* Vertex */
56957459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
57057459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
57157459e9aSMatthew G Knepley     /* X Face */
57282f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
573201c01f8SBarry Smith     else if (dim == 2) {
57457459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
57557459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
576da80777bSKarl Rupp       PetscInt points[3];
57757459e9aSMatthew G Knepley 
578da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
57957459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
58082f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
58157459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
58257459e9aSMatthew G Knepley     /* Y Face */
58382f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
584201c01f8SBarry Smith     else if (dim == 2) {
58557459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
58657459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
587da80777bSKarl Rupp       PetscInt points[3];
58857459e9aSMatthew G Knepley 
589da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
59057459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
59182f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
59257459e9aSMatthew G Knepley   } else {
59357459e9aSMatthew G Knepley     /* Z Face */
59482f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
59582f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
59682f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
59757459e9aSMatthew G Knepley   }
598aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
599aa1993deSMatthew G Knepley }
600aa1993deSMatthew G Knepley 
601aa1993deSMatthew G Knepley #undef __FUNCT__
602aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
603aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
604aa1993deSMatthew G Knepley {
605aa1993deSMatthew G Knepley   PetscScalar    *vArray;
606aa1993deSMatthew G Knepley   PetscErrorCode ierr;
607aa1993deSMatthew G Knepley 
608aa1993deSMatthew G Knepley   PetscFunctionBegin;
609aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
610aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
611aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
612aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
613aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
61457459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
61557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
61657459e9aSMatthew G Knepley }
61757459e9aSMatthew G Knepley 
618ed4e918cSMatthew G Knepley #undef __FUNCT__
619ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
620ed4e918cSMatthew G Knepley /*@
621f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
622341c9bdaSMatthew G Knepley 
623ed4e918cSMatthew G Knepley   Not Collective
624ed4e918cSMatthew G Knepley 
625ed4e918cSMatthew G Knepley   Input Parameter:
626ed4e918cSMatthew G Knepley + da - the distributed array
627ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
628ed4e918cSMatthew G Knepley 
629ed4e918cSMatthew G Knepley   Output Parameter:
630ed4e918cSMatthew G Knepley . cell - the local cell number
631ed4e918cSMatthew G Knepley 
632ed4e918cSMatthew G Knepley   Level: developer
633ed4e918cSMatthew G Knepley 
634ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
635ed4e918cSMatthew G Knepley @*/
636ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
637341c9bdaSMatthew G Knepley {
638341c9bdaSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
639*c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
6404e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
641ed4e918cSMatthew 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;
642341c9bdaSMatthew G Knepley 
643341c9bdaSMatthew G Knepley   PetscFunctionBegin;
644ed4e918cSMatthew G Knepley   *cell = -1;
64582f516ccSBarry Smith   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w);
64682f516ccSBarry Smith   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye);
64782f516ccSBarry Smith   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze);
648ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
649ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
650341c9bdaSMatthew G Knepley }
651341c9bdaSMatthew G Knepley 
65257459e9aSMatthew G Knepley #undef __FUNCT__
65357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
65415d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
65557459e9aSMatthew G Knepley {
65657459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
65757459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
65857459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
65957459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
66057459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
66157459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
66257459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
66357459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
66457459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
66557459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
66657459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
66757459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
66857459e9aSMatthew G Knepley   PetscReal         invDet;
66957459e9aSMatthew G Knepley   PetscErrorCode    ierr;
67057459e9aSMatthew G Knepley 
67157459e9aSMatthew G Knepley   PetscFunctionBegin;
672e477d84eSMatthew G. Knepley #if 0
67315d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
6744c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
67515d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
67615d2e57cSJed Brown #endif
67715d2e57cSJed Brown   J[0]    = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
67815d2e57cSJed Brown   J[2]    = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
67957459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
68057459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
681e477d84eSMatthew G. Knepley   if (invJ) {
68257459e9aSMatthew G Knepley     invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
68357459e9aSMatthew G Knepley     invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
684e477d84eSMatthew G. Knepley   }
68557459e9aSMatthew G Knepley   ierr    = PetscLogFlops(30);CHKERRQ(ierr);
68657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
68757459e9aSMatthew G Knepley }
68857459e9aSMatthew G Knepley 
68957459e9aSMatthew G Knepley #undef __FUNCT__
69057459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
69121454ff5SMatthew G. Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
69257459e9aSMatthew G Knepley {
69357459e9aSMatthew G Knepley   DM               cdm;
69457459e9aSMatthew G Knepley   Vec              coordinates;
69521454ff5SMatthew G. Knepley   const PetscReal *quadPoints;
696d7a12f93SMatthew G. Knepley   PetscScalar     *vertices = NULL;
69721454ff5SMatthew G. Knepley   PetscInt         numQuadPoints, csize, dim, d, q;
69857459e9aSMatthew G Knepley   PetscErrorCode   ierr;
69957459e9aSMatthew G Knepley 
70057459e9aSMatthew G Knepley   PetscFunctionBegin;
70157459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
70257459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
703e477d84eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
7046636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
7050a3ada39SMatthew G. Knepley   ierr = DMDAVecGetClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
7068865f1eaSKarl Rupp   for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]);
70757459e9aSMatthew G Knepley   switch (dim) {
70857459e9aSMatthew G Knepley   case 2:
70921454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, NULL);CHKERRQ(ierr);
71021454ff5SMatthew G. Knepley     for (q = 0; q < numQuadPoints; ++q) {
71121454ff5SMatthew G. Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
71257459e9aSMatthew G Knepley     }
71357459e9aSMatthew G Knepley     break;
71457459e9aSMatthew G Knepley   default:
71582f516ccSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension %d not supported", dim);
71657459e9aSMatthew G Knepley   }
7170a3ada39SMatthew G. Knepley   ierr = DMDAVecRestoreClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
71857459e9aSMatthew G Knepley   PetscFunctionReturn(0);
71957459e9aSMatthew G Knepley }
720