xref: /petsc/src/dm/field/impls/ds/dmfieldds.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
20145028aSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/
351a454edSToby Isaac #include <petscfe.h>
451a454edSToby Isaac #include <petscdmplex.h>
551a454edSToby Isaac #include <petscds.h>
651a454edSToby Isaac 
751a454edSToby Isaac typedef struct _n_DMField_DS
851a454edSToby Isaac {
951a454edSToby Isaac   PetscInt    fieldNum;
1051a454edSToby Isaac   Vec         vec;
1151a454edSToby Isaac   PetscInt    height;
1251a454edSToby Isaac   PetscObject *disc;
1351a454edSToby Isaac   PetscBool   multifieldVec;
1451a454edSToby Isaac }
1551a454edSToby Isaac DMField_DS;
1651a454edSToby Isaac 
1751a454edSToby Isaac static PetscErrorCode DMFieldDestroy_DS(DMField field)
1851a454edSToby Isaac {
1951a454edSToby Isaac   DMField_DS     *dsfield;
2051a454edSToby Isaac   PetscInt       i;
2151a454edSToby Isaac 
2251a454edSToby Isaac   PetscFunctionBegin;
2351a454edSToby Isaac   dsfield = (DMField_DS *) field->data;
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&dsfield->vec));
2551a454edSToby Isaac   for (i = 0; i < dsfield->height; i++) {
26*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectDereference(dsfield->disc[i]));
2751a454edSToby Isaac   }
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(dsfield->disc));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(dsfield));
3051a454edSToby Isaac   PetscFunctionReturn(0);
3151a454edSToby Isaac }
3251a454edSToby Isaac 
3351a454edSToby Isaac static PetscErrorCode DMFieldView_DS(DMField field,PetscViewer viewer)
3451a454edSToby Isaac {
3551a454edSToby Isaac   DMField_DS     *dsfield = (DMField_DS *) field->data;
3651a454edSToby Isaac   PetscBool      iascii;
3751a454edSToby Isaac   PetscObject    disc;
3851a454edSToby Isaac 
3951a454edSToby Isaac   PetscFunctionBegin;
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
4151a454edSToby Isaac   disc = dsfield->disc[0];
4251a454edSToby Isaac   if (iascii) {
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer, "PetscDS field %D\n",dsfield->fieldNum));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectView(disc,viewer));
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopTab(viewer));
4751a454edSToby Isaac   }
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPushTab(viewer));
492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dsfield->multifieldVec,PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"View of subfield not implemented yet");
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(dsfield->vec,viewer));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPopTab(viewer));
5251a454edSToby Isaac   PetscFunctionReturn(0);
5351a454edSToby Isaac }
5451a454edSToby Isaac 
5551a454edSToby Isaac static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject *disc)
5651a454edSToby Isaac {
5751a454edSToby Isaac   DMField_DS     *dsfield = (DMField_DS *) field->data;
5851a454edSToby Isaac 
5951a454edSToby Isaac   PetscFunctionBegin;
6051a454edSToby Isaac   if (!dsfield->disc[height]) {
6151a454edSToby Isaac     PetscClassId   id;
6251a454edSToby Isaac 
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(dsfield->disc[0],&id));
6451a454edSToby Isaac     if (id == PETSCFE_CLASSID) {
6551a454edSToby Isaac       PetscFE fe = (PetscFE) dsfield->disc[0];
6651a454edSToby Isaac 
67*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFECreateHeightTrace(fe,height,(PetscFE *)&dsfield->disc[height]));
6851a454edSToby Isaac     }
6951a454edSToby Isaac   }
7051a454edSToby Isaac   *disc = dsfield->disc[height];
7151a454edSToby Isaac   PetscFunctionReturn(0);
7251a454edSToby Isaac }
7351a454edSToby Isaac 
7462bd480fSMatthew G. Knepley /* y[m,c] = A[m,n,c] . b[n] */
7551a454edSToby Isaac #define DMFieldDSdot(y,A,b,m,n,c,cast)                                           \
7651a454edSToby Isaac   do {                                                                           \
7751a454edSToby Isaac     PetscInt _i, _j, _k;                                                         \
7851a454edSToby Isaac     for (_i = 0; _i < (m); _i++) {                                               \
7951a454edSToby Isaac       for (_k = 0; _k < (c); _k++) {                                             \
8051a454edSToby Isaac         (y)[_i * (c) + _k] = 0.;                                                 \
8151a454edSToby Isaac       }                                                                          \
8251a454edSToby Isaac       for (_j = 0; _j < (n); _j++) {                                             \
8351a454edSToby Isaac         for (_k = 0; _k < (c); _k++) {                                           \
8451a454edSToby Isaac           (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \
8551a454edSToby Isaac         }                                                                        \
8651a454edSToby Isaac       }                                                                          \
8751a454edSToby Isaac     }                                                                            \
8851a454edSToby Isaac   } while (0)
8951a454edSToby Isaac 
90ef0bb6c7SMatthew G. Knepley /* TODO: Reorganize interface so that I can reuse a tabulation rather than mallocing each time */
9151a454edSToby Isaac static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H)
9251a454edSToby Isaac {
9351a454edSToby Isaac   DMField_DS      *dsfield = (DMField_DS *) field->data;
9451a454edSToby Isaac   DM              dm;
9551a454edSToby Isaac   PetscObject     disc;
9651a454edSToby Isaac   PetscClassId    classid;
9751a454edSToby Isaac   PetscInt        nq, nc, dim, meshDim, numCells;
9851a454edSToby Isaac   PetscSection    section;
9951a454edSToby Isaac   const PetscReal *qpoints;
10051a454edSToby Isaac   PetscBool       isStride;
10151a454edSToby Isaac   const PetscInt  *points = NULL;
10251a454edSToby Isaac   PetscInt        sfirst = -1, stride = -1;
10351a454edSToby Isaac 
10451a454edSToby Isaac   PetscFunctionBeginHot;
10551a454edSToby Isaac   dm   = field->dm;
10651a454edSToby Isaac   nc   = field->numComponents;
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(quad,&dim,NULL,&nq,&qpoints,NULL));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldDSGetHeightDisc(field,dsfield->height - 1 - dim,&disc));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&meshDim));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm,&section));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetField(section,dsfield->fieldNum,&section));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetClassId(disc,&classid));
11351a454edSToby Isaac   /* TODO: batch */
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pointIS,ISSTRIDE,&isStride));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(pointIS,&numCells));
11651a454edSToby Isaac   if (isStride) {
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISStrideGetInfo(pointIS,&sfirst,&stride));
11851a454edSToby Isaac   } else {
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(pointIS,&points));
12051a454edSToby Isaac   }
12151a454edSToby Isaac   if (classid == PETSCFE_CLASSID) {
12251a454edSToby Isaac     PetscFE      fe = (PetscFE) disc;
12351a454edSToby Isaac     PetscInt     feDim, i;
124ef0bb6c7SMatthew G. Knepley     PetscInt          K = H ? 2 : (D ? 1 : (B ? 0 : -1));
125ef0bb6c7SMatthew G. Knepley     PetscTabulation   T;
12651a454edSToby Isaac 
127*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetDimension(fe,&feDim));
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECreateTabulation(fe,1,nq,qpoints,K,&T));
12951a454edSToby Isaac     for (i = 0; i < numCells; i++) {
13051a454edSToby Isaac       PetscInt     c = isStride ? (sfirst + i * stride) : points[i];
1312710589cSToby Isaac       PetscInt     closureSize;
1322710589cSToby Isaac       PetscScalar *elem = NULL;
13351a454edSToby Isaac 
134*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecGetClosure(dm,section,dsfield->vec,c,&closureSize,&elem));
13551a454edSToby Isaac       if (B) {
13662bd480fSMatthew G. Knepley         /* field[c] = T[q,b,c] . coef[b], so v[c] = T[q,b,c] . coords[b] */
13751a454edSToby Isaac         if (type == PETSC_SCALAR) {
13851a454edSToby Isaac           PetscScalar *cB = &((PetscScalar *) B)[nc * nq * i];
13951a454edSToby Isaac 
140ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,(PetscScalar));
14151a454edSToby Isaac         } else {
14251a454edSToby Isaac           PetscReal *cB = &((PetscReal *) B)[nc * nq * i];
14351a454edSToby Isaac 
144ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,PetscRealPart);
14551a454edSToby Isaac         }
14651a454edSToby Isaac       }
14751a454edSToby Isaac       if (D) {
14851a454edSToby Isaac         if (type == PETSC_SCALAR) {
14951a454edSToby Isaac           PetscScalar *cD = &((PetscScalar *) D)[nc * nq * dim * i];
15051a454edSToby Isaac 
151ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),(PetscScalar));
15251a454edSToby Isaac         } else {
15351a454edSToby Isaac           PetscReal *cD = &((PetscReal *) D)[nc * nq * dim * i];
15451a454edSToby Isaac 
155ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),PetscRealPart);
15651a454edSToby Isaac         }
15751a454edSToby Isaac       }
15851a454edSToby Isaac       if (H) {
15951a454edSToby Isaac         if (type == PETSC_SCALAR) {
16051a454edSToby Isaac           PetscScalar *cH = &((PetscScalar *) H)[nc * nq * dim * dim * i];
16151a454edSToby Isaac 
162ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),(PetscScalar));
16351a454edSToby Isaac         } else {
16451a454edSToby Isaac           PetscReal *cH = &((PetscReal *) H)[nc * nq * dim * dim * i];
16551a454edSToby Isaac 
166ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),PetscRealPart);
16751a454edSToby Isaac         }
16851a454edSToby Isaac       }
169*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecRestoreClosure(dm,section,dsfield->vec,c,&closureSize,&elem));
17051a454edSToby Isaac     }
171*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscTabulationDestroy(&T));
1722da392ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented");
17351a454edSToby Isaac   if (!isStride) {
174*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(pointIS,&points));
17551a454edSToby Isaac   }
17651a454edSToby Isaac   PetscFunctionReturn(0);
17751a454edSToby Isaac }
17851a454edSToby Isaac 
179a6216909SToby Isaac static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
180a6216909SToby Isaac {
181a6216909SToby Isaac   DMField_DS        *dsfield = (DMField_DS *) field->data;
182a6216909SToby Isaac   PetscSF            cellSF = NULL;
183a6216909SToby Isaac   const PetscSFNode *cells;
184a6216909SToby Isaac   PetscInt           c, nFound, numCells, feDim, nc;
185a6216909SToby Isaac   const PetscInt    *cellDegrees;
186a6216909SToby Isaac   const PetscScalar *pointsArray;
187a6216909SToby Isaac   PetscScalar       *cellPoints;
188a6216909SToby Isaac   PetscInt           gatherSize, gatherMax;
189a6216909SToby Isaac   PetscInt           dim, dimR, offset;
190a6216909SToby Isaac   MPI_Datatype       pointType;
191a6216909SToby Isaac   PetscObject        cellDisc;
192a6216909SToby Isaac   PetscFE            cellFE;
193a6216909SToby Isaac   PetscClassId       discID;
194a6216909SToby Isaac   PetscReal         *coordsReal, *coordsRef;
195a6216909SToby Isaac   PetscSection       section;
196a6216909SToby Isaac   PetscScalar       *cellBs = NULL, *cellDs = NULL, *cellHs = NULL;
197a6216909SToby Isaac   PetscReal         *cellBr = NULL, *cellDr = NULL, *cellHr = NULL;
1984cbe5319SToby Isaac   PetscReal         *v, *J, *invJ, *detJ;
199a6216909SToby Isaac 
200a6216909SToby Isaac   PetscFunctionBegin;
201a6216909SToby Isaac   nc   = field->numComponents;
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(field->dm,&section));
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldDSGetHeightDisc(field,0,&cellDisc));
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetClassId(cellDisc, &discID));
2052c71b3e2SJacob Faibussowitsch   PetscCheckFalse(discID != PETSCFE_CLASSID,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Discretization type not supported");
206a6216909SToby Isaac   cellFE = (PetscFE) cellDisc;
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDimension(cellFE,&feDim));
208*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(field->dm, &dim));
209*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(field->dm, &dimR));
210*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF));
211*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells));
212a6216909SToby Isaac   for (c = 0; c < nFound; c++) {
2132c71b3e2SJacob Faibussowitsch     PetscCheckFalse(cells[c].index < 0,PetscObjectComm((PetscObject)points),PETSC_ERR_ARG_WRONG, "Point %D could not be located", c);
214a6216909SToby Isaac   }
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFComputeDegreeBegin(cellSF,&cellDegrees));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFComputeDegreeEnd(cellSF,&cellDegrees));
217a6216909SToby Isaac   for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) {
218a6216909SToby Isaac     gatherMax = PetscMax(gatherMax,cellDegrees[c]);
219a6216909SToby Isaac     gatherSize += cellDegrees[c];
220a6216909SToby Isaac   }
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(gatherSize*dim,&cellPoints,gatherMax*dim,&coordsReal,gatherMax*dimR,&coordsRef));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(gatherMax*dimR,&v,gatherMax*dimR*dimR,&J,gatherMax*dimR*dimR,&invJ,gatherMax,&detJ));
223a6216909SToby Isaac   if (datatype == PETSC_SCALAR) {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs));
225a6216909SToby Isaac   } else {
226*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr));
227a6216909SToby Isaac   }
228a6216909SToby Isaac 
229*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_contiguous(dim,MPIU_SCALAR,&pointType));
230*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_commit(&pointType));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(points,&pointsArray));
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints));
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(points,&pointsArray));
235a6216909SToby Isaac   for (c = 0, offset = 0; c < numCells; c++) {
236a6216909SToby Isaac     PetscInt nq = cellDegrees[c], p;
237a6216909SToby Isaac 
238a6216909SToby Isaac     if (nq) {
239ef0bb6c7SMatthew G. Knepley       PetscInt          K = H ? 2 : (D ? 1 : (B ? 0 : -1));
240ef0bb6c7SMatthew G. Knepley       PetscTabulation   T;
241a6216909SToby Isaac       PetscInt     closureSize;
242a6216909SToby Isaac       PetscScalar *elem = NULL;
2434cbe5319SToby Isaac       PetscReal   *quadPoints;
2444cbe5319SToby Isaac       PetscQuadrature quad;
2454cbe5319SToby Isaac       PetscInt d, e, f, g;
246a6216909SToby Isaac 
2474d1a973fSToby Isaac       for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]);
248*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef));
249*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFECreateTabulation(cellFE,1,nq,coordsRef,K,&T));
250*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureCreate(PETSC_COMM_SELF, &quad));
251*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(dimR * nq, &quadPoints));
2524cbe5319SToby Isaac       for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p];
253*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL));
254*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ));
255*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureDestroy(&quad));
256*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecGetClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem));
257a6216909SToby Isaac       if (B) {
258a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
259a6216909SToby Isaac           PetscScalar *cB = &cellBs[nc * offset];
260a6216909SToby Isaac 
261ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,(PetscScalar));
262a6216909SToby Isaac         } else {
263a6216909SToby Isaac           PetscReal *cB = &cellBr[nc * offset];
264a6216909SToby Isaac 
265ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,PetscRealPart);
266a6216909SToby Isaac         }
267a6216909SToby Isaac       }
268a6216909SToby Isaac       if (D) {
269a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
270a6216909SToby Isaac           PetscScalar *cD = &cellDs[nc * dim * offset];
271a6216909SToby Isaac 
272ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),(PetscScalar));
2734cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
2744cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
2754d1a973fSToby Isaac               PetscScalar vs[3];
2764d1a973fSToby Isaac 
2774cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
2784d1a973fSToby Isaac                 vs[d] = 0.;
2794cbe5319SToby Isaac                 for (e = 0; e < dimR; e++) {
2804d1a973fSToby Isaac                   vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
2814cbe5319SToby Isaac                 }
2824cbe5319SToby Isaac               }
2834cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
2844d1a973fSToby Isaac                 cD[(nc * p + g) * dimR + d] = vs[d];
2854cbe5319SToby Isaac               }
2864cbe5319SToby Isaac             }
2874cbe5319SToby Isaac           }
288a6216909SToby Isaac         } else {
289a6216909SToby Isaac           PetscReal *cD = &cellDr[nc * dim * offset];
290a6216909SToby Isaac 
291ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),PetscRealPart);
2924cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
2934cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
2944cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
2954cbe5319SToby Isaac                 v[d] = 0.;
2964cbe5319SToby Isaac                 for (e = 0; e < dimR; e++) {
2974cbe5319SToby Isaac                   v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
2984cbe5319SToby Isaac                 }
2994cbe5319SToby Isaac               }
3004cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3014cbe5319SToby Isaac                 cD[(nc * p + g) * dimR + d] = v[d];
3024cbe5319SToby Isaac               }
3034cbe5319SToby Isaac             }
3044cbe5319SToby Isaac           }
305a6216909SToby Isaac         }
306a6216909SToby Isaac       }
307a6216909SToby Isaac       if (H) {
308a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
309a6216909SToby Isaac           PetscScalar *cH = &cellHs[nc * dim * dim * offset];
310a6216909SToby Isaac 
311ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),(PetscScalar));
3124cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3134cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3144d1a973fSToby Isaac               PetscScalar vs[3];
3154d1a973fSToby Isaac 
3164cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3174d1a973fSToby Isaac                 vs[d] = 0.;
3184cbe5319SToby Isaac                 for (e = 0; e < dimR; e++) {
3194d1a973fSToby Isaac                   vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3204cbe5319SToby Isaac                 }
3214cbe5319SToby Isaac               }
3224cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3234d1a973fSToby Isaac                 cH[(nc * dimR * p + g) * dimR + d] = vs[d];
3244cbe5319SToby Isaac               }
3254cbe5319SToby Isaac             }
3264cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3274cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3284d1a973fSToby Isaac                 PetscScalar vs[3];
3294d1a973fSToby Isaac 
3304cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3314d1a973fSToby Isaac                   vs[d] = 0.;
3324cbe5319SToby Isaac                   for (e = 0; e < dimR; e++) {
3334d1a973fSToby Isaac                     vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3344cbe5319SToby Isaac                   }
3354cbe5319SToby Isaac                 }
3364cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3374d1a973fSToby Isaac                   cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d];
3384cbe5319SToby Isaac                 }
3394cbe5319SToby Isaac               }
3404cbe5319SToby Isaac             }
3414cbe5319SToby Isaac           }
342a6216909SToby Isaac         } else {
343a6216909SToby Isaac           PetscReal *cH = &cellHr[nc * dim * dim * offset];
344a6216909SToby Isaac 
345ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),PetscRealPart);
3464cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3474cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3484cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3494cbe5319SToby Isaac                 v[d] = 0.;
3504cbe5319SToby Isaac                 for (e = 0; e < dimR; e++) {
3514cbe5319SToby Isaac                   v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3524cbe5319SToby Isaac                 }
3534cbe5319SToby Isaac               }
3544cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3554cbe5319SToby Isaac                 cH[(nc * dimR * p + g) * dimR + d] = v[d];
3564cbe5319SToby Isaac               }
3574cbe5319SToby Isaac             }
3584cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3594cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3604cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3614cbe5319SToby Isaac                   v[d] = 0.;
3624cbe5319SToby Isaac                   for (e = 0; e < dimR; e++) {
3634cbe5319SToby Isaac                     v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3644cbe5319SToby Isaac                   }
3654cbe5319SToby Isaac                 }
3664cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3674cbe5319SToby Isaac                   cH[((nc * p + g) * dimR + d) * dimR + f] = v[d];
3684cbe5319SToby Isaac                 }
3694cbe5319SToby Isaac               }
3704cbe5319SToby Isaac             }
3714cbe5319SToby Isaac           }
372a6216909SToby Isaac         }
373a6216909SToby Isaac       }
374*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecRestoreClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem));
375*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscTabulationDestroy(&T));
376a6216909SToby Isaac     }
377a6216909SToby Isaac     offset += nq;
378a6216909SToby Isaac   }
379a6216909SToby Isaac   {
380a6216909SToby Isaac     MPI_Datatype origtype;
3814cbe5319SToby Isaac 
382a6216909SToby Isaac     if (datatype == PETSC_SCALAR) {
383a6216909SToby Isaac       origtype = MPIU_SCALAR;
384a6216909SToby Isaac     } else {
385a6216909SToby Isaac       origtype = MPIU_REAL;
386a6216909SToby Isaac     }
387a6216909SToby Isaac     if (B) {
388a6216909SToby Isaac       MPI_Datatype Btype;
389a6216909SToby Isaac 
390*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_contiguous(nc, origtype, &Btype));
391*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_commit(&Btype));
392*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterBegin(cellSF,Btype,(datatype == PETSC_SCALAR) ? (void *) cellBs : (void *) cellBr, B));
393*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterEnd(cellSF,Btype,(datatype == PETSC_SCALAR) ? (void *) cellBs : (void *) cellBr, B));
394*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_free(&Btype));
395a6216909SToby Isaac     }
396a6216909SToby Isaac     if (D) {
397a6216909SToby Isaac       MPI_Datatype Dtype;
398a6216909SToby Isaac 
399*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_contiguous(nc * dim, origtype, &Dtype));
400*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_commit(&Dtype));
401*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterBegin(cellSF,Dtype,(datatype == PETSC_SCALAR) ? (void *) cellDs : (void *) cellDr, D));
402*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterEnd(cellSF,Dtype,(datatype == PETSC_SCALAR) ? (void *) cellDs : (void *) cellDr, D));
403*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_free(&Dtype));
404a6216909SToby Isaac     }
405a6216909SToby Isaac     if (H) {
406a6216909SToby Isaac       MPI_Datatype Htype;
407a6216909SToby Isaac 
408*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_contiguous(nc * dim * dim, origtype, &Htype));
409*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_commit(&Htype));
410*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterBegin(cellSF,Htype,(datatype == PETSC_SCALAR) ? (void *) cellHs : (void *) cellHr, H));
411*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFScatterEnd(cellSF,Htype,(datatype == PETSC_SCALAR) ? (void *) cellHs : (void *) cellHr, H));
412*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Type_free(&Htype));
413a6216909SToby Isaac     }
414a6216909SToby Isaac   }
415*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(v,J,invJ,detJ));
416*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(cellBr, cellDr, cellHr));
417*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(cellBs, cellDs, cellHs));
418*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(cellPoints,coordsReal,coordsRef));
419*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_free(&pointType));
420*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&cellSF));
421a6216909SToby Isaac   PetscFunctionReturn(0);
422a6216909SToby Isaac }
423a6216909SToby Isaac 
424805e7170SToby Isaac static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H)
425805e7170SToby Isaac {
426805e7170SToby Isaac   DMField_DS      *dsfield = (DMField_DS *) field->data;
427805e7170SToby Isaac   PetscInt         h, imin;
428805e7170SToby Isaac   PetscInt         dim;
429805e7170SToby Isaac   PetscClassId     id;
430805e7170SToby Isaac   PetscQuadrature  quad = NULL;
431b7260050SToby Isaac   PetscInt         maxDegree;
432805e7170SToby Isaac   PetscFEGeom      *geom;
433805e7170SToby Isaac   PetscInt         Nq, Nc, dimC, qNc, N;
434805e7170SToby Isaac   PetscInt         numPoints;
435805e7170SToby Isaac   void            *qB = NULL, *qD = NULL, *qH = NULL;
436805e7170SToby Isaac   const PetscReal *weights;
437805e7170SToby Isaac   MPI_Datatype     mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL;
438805e7170SToby Isaac   PetscObject      disc;
439805e7170SToby Isaac   DMField          coordField;
440805e7170SToby Isaac 
441805e7170SToby Isaac   PetscFunctionBegin;
442805e7170SToby Isaac   Nc = field->numComponents;
443*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(field->dm, &dimC));
444*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(field->dm, &dim));
445*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(pointIS, &numPoints));
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetMinMax(pointIS,&imin,NULL));
447805e7170SToby Isaac   for (h = 0; h < dsfield->height; h++) {
448805e7170SToby Isaac     PetscInt hEnd;
449805e7170SToby Isaac 
450*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd));
451805e7170SToby Isaac     if (imin < hEnd) break;
452805e7170SToby Isaac   }
453805e7170SToby Isaac   dim -= h;
454*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldDSGetHeightDisc(field,h,&disc));
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetClassId(disc,&id));
4562c71b3e2SJacob Faibussowitsch   PetscCheckFalse(id != PETSCFE_CLASSID,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported");
457*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateField(field->dm, &coordField));
458*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
459b7260050SToby Isaac   if (maxDegree <= 1) {
460*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
461805e7170SToby Isaac   }
462*5f80ce2aSJacob Faibussowitsch   if (!quad) CHKERRQ(DMFieldCreateDefaultQuadrature(field, pointIS, &quad));
4632c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!quad,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages");
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&geom));
465*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights));
4662c71b3e2SJacob Faibussowitsch   PetscCheckFalse(qNc != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components");
467805e7170SToby Isaac   N = numPoints * Nq * Nc;
468*5f80ce2aSJacob Faibussowitsch   if (B) CHKERRQ(DMGetWorkArray(field->dm, N, mpitype, &qB));
469*5f80ce2aSJacob Faibussowitsch   if (D) CHKERRQ(DMGetWorkArray(field->dm, N * dimC, mpitype, &qD));
470*5f80ce2aSJacob Faibussowitsch   if (H) CHKERRQ(DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
471*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldEvaluateFE(field,pointIS,quad,type,qB,qD,qH));
472805e7170SToby Isaac   if (B) {
473805e7170SToby Isaac     PetscInt i, j, k;
474805e7170SToby Isaac 
475805e7170SToby Isaac     if (type == PETSC_SCALAR) {
476805e7170SToby Isaac       PetscScalar * sB  = (PetscScalar *) B;
477805e7170SToby Isaac       PetscScalar * sqB = (PetscScalar *) qB;
478805e7170SToby Isaac 
479805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
480805e7170SToby Isaac         PetscReal vol = 0.;
481805e7170SToby Isaac 
482805e7170SToby Isaac         for (j = 0; j < Nc; j++) {sB[i * Nc + j] = 0.;}
483805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
484805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
485805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
486805e7170SToby Isaac             sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[ (i * Nq + k) * Nc + j];
487805e7170SToby Isaac           }
488805e7170SToby Isaac         }
48962bd480fSMatthew G. Knepley         for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol;
490805e7170SToby Isaac       }
491805e7170SToby Isaac     } else {
492805e7170SToby Isaac       PetscReal * rB  = (PetscReal *) B;
493805e7170SToby Isaac       PetscReal * rqB = (PetscReal *) qB;
494805e7170SToby Isaac 
495805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
496805e7170SToby Isaac         PetscReal vol = 0.;
497805e7170SToby Isaac 
498805e7170SToby Isaac         for (j = 0; j < Nc; j++) {rB[i * Nc + j] = 0.;}
499805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
500805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
501805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
502805e7170SToby Isaac             rB[i * Nc + j] += weights[k] * rqB[ (i * Nq + k) * Nc + j];
503805e7170SToby Isaac           }
504805e7170SToby Isaac         }
505805e7170SToby Isaac         for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol;
506805e7170SToby Isaac       }
507805e7170SToby Isaac     }
508805e7170SToby Isaac   }
509805e7170SToby Isaac   if (D) {
510805e7170SToby Isaac     PetscInt i, j, k, l, m;
511805e7170SToby Isaac 
512805e7170SToby Isaac     if (type == PETSC_SCALAR) {
513805e7170SToby Isaac       PetscScalar * sD  = (PetscScalar *) D;
514805e7170SToby Isaac       PetscScalar * sqD = (PetscScalar *) qD;
515805e7170SToby Isaac 
516805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
517805e7170SToby Isaac         PetscReal vol = 0.;
518805e7170SToby Isaac 
519805e7170SToby Isaac         for (j = 0; j < Nc * dimC; j++) {sD[i * Nc * dimC + j] = 0.;}
520805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
521805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
522805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
523805e7170SToby Isaac             PetscScalar pD[3] = {0.,0.,0.};
524805e7170SToby Isaac 
525805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
526805e7170SToby Isaac               for (m = 0; m < dim; m++) {
527805e7170SToby Isaac                 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * sqD[((i * Nq + k) * Nc + j) * dim + m];
528805e7170SToby Isaac               }
529805e7170SToby Isaac             }
530805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
531805e7170SToby Isaac               sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
532805e7170SToby Isaac             }
533805e7170SToby Isaac           }
534805e7170SToby Isaac         }
535805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol;
536805e7170SToby Isaac       }
537805e7170SToby Isaac     } else {
538805e7170SToby Isaac       PetscReal * rD  = (PetscReal *) D;
539805e7170SToby Isaac       PetscReal * rqD = (PetscReal *) qD;
540805e7170SToby Isaac 
541805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
542805e7170SToby Isaac         PetscReal vol = 0.;
543805e7170SToby Isaac 
544805e7170SToby Isaac         for (j = 0; j < Nc * dimC; j++) {rD[i * Nc * dimC + j] = 0.;}
545805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
546805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
547805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
548805e7170SToby Isaac             PetscReal pD[3] = {0.,0.,0.};
549805e7170SToby Isaac 
550805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
551805e7170SToby Isaac               for (m = 0; m < dim; m++) {
552805e7170SToby Isaac                 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * rqD[((i * Nq + k) * Nc + j) * dim + m];
553805e7170SToby Isaac               }
554805e7170SToby Isaac             }
555805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
556805e7170SToby Isaac               rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
557805e7170SToby Isaac             }
558805e7170SToby Isaac           }
559805e7170SToby Isaac         }
560805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol;
561805e7170SToby Isaac       }
562805e7170SToby Isaac     }
563805e7170SToby Isaac   }
564805e7170SToby Isaac   if (H) {
565805e7170SToby Isaac     PetscInt i, j, k, l, m, q, r;
566805e7170SToby Isaac 
567805e7170SToby Isaac     if (type == PETSC_SCALAR) {
568805e7170SToby Isaac       PetscScalar * sH  = (PetscScalar *) H;
569805e7170SToby Isaac       PetscScalar * sqH = (PetscScalar *) qH;
570805e7170SToby Isaac 
571805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
572805e7170SToby Isaac         PetscReal vol = 0.;
573805e7170SToby Isaac 
574805e7170SToby Isaac         for (j = 0; j < Nc * dimC * dimC; j++) {sH[i * Nc * dimC * dimC + j] = 0.;}
575805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
5764d1a973fSToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
577805e7170SToby Isaac 
578805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
579805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
580805e7170SToby Isaac             PetscScalar pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}};
581805e7170SToby Isaac             const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC];
582805e7170SToby Isaac 
583805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
584805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
585805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
586805e7170SToby Isaac                   for (r = 0; r < dim; r++) {
587805e7170SToby Isaac                     pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r];
588805e7170SToby Isaac                   }
589805e7170SToby Isaac                 }
590805e7170SToby Isaac               }
591805e7170SToby Isaac             }
592805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
593805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
594805e7170SToby Isaac                 sH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
595805e7170SToby Isaac               }
596805e7170SToby Isaac             }
597805e7170SToby Isaac           }
598805e7170SToby Isaac         }
599805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol;
600805e7170SToby Isaac       }
601805e7170SToby Isaac     } else {
602805e7170SToby Isaac       PetscReal * rH  = (PetscReal *) H;
603805e7170SToby Isaac       PetscReal * rqH = (PetscReal *) qH;
604805e7170SToby Isaac 
605805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
606805e7170SToby Isaac         PetscReal vol = 0.;
607805e7170SToby Isaac 
608805e7170SToby Isaac         for (j = 0; j < Nc * dimC * dimC; j++) {rH[i * Nc * dimC * dimC + j] = 0.;}
609805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
610805e7170SToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
611805e7170SToby Isaac 
612805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
613805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
614805e7170SToby Isaac             PetscReal pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}};
615805e7170SToby Isaac             const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC];
616805e7170SToby Isaac 
617805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
618805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
619805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
620805e7170SToby Isaac                   for (r = 0; r < dim; r++) {
621805e7170SToby Isaac                     pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r];
622805e7170SToby Isaac                   }
623805e7170SToby Isaac                 }
624805e7170SToby Isaac               }
625805e7170SToby Isaac             }
626805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
627805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
628805e7170SToby Isaac                 rH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
629805e7170SToby Isaac               }
630805e7170SToby Isaac             }
631805e7170SToby Isaac           }
632805e7170SToby Isaac         }
633805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol;
634805e7170SToby Isaac       }
635805e7170SToby Isaac     }
636805e7170SToby Isaac   }
637*5f80ce2aSJacob Faibussowitsch   if (B) CHKERRQ(DMRestoreWorkArray(field->dm, N, mpitype, &qB));
638*5f80ce2aSJacob Faibussowitsch   if (D) CHKERRQ(DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD));
639*5f80ce2aSJacob Faibussowitsch   if (H) CHKERRQ(DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
640*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGeomDestroy(&geom));
641*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&quad));
642805e7170SToby Isaac   PetscFunctionReturn(0);
643805e7170SToby Isaac }
644805e7170SToby Isaac 
645b7260050SToby Isaac static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
64651a454edSToby Isaac {
64751a454edSToby Isaac   DMField_DS     *dsfield;
64851a454edSToby Isaac   PetscObject    disc;
649741db25dSToby Isaac   PetscInt       h, imin, imax;
65051a454edSToby Isaac   PetscClassId   id;
65151a454edSToby Isaac 
65251a454edSToby Isaac   PetscFunctionBegin;
65351a454edSToby Isaac   dsfield = (DMField_DS *) field->data;
654*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetMinMax(pointIS,&imin,&imax));
655741db25dSToby Isaac   if (imin >= imax) {
65651aa0bf7SToby Isaac     h = 0;
65751aa0bf7SToby Isaac   } else {
658f99c8401SToby Isaac     for (h = 0; h < dsfield->height; h++) {
65951a454edSToby Isaac       PetscInt hEnd;
66051a454edSToby Isaac 
661*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd));
66251a454edSToby Isaac       if (imin < hEnd) break;
66351a454edSToby Isaac     }
66451aa0bf7SToby Isaac   }
665*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldDSGetHeightDisc(field,h,&disc));
666*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetClassId(disc,&id));
66751a454edSToby Isaac   if (id == PETSCFE_CLASSID) {
66851a454edSToby Isaac     PetscFE    fe = (PetscFE) disc;
66951a454edSToby Isaac     PetscSpace sp;
67051a454edSToby Isaac 
671*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetBasisSpace(fe, &sp));
672*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSpaceGetDegree(sp, minDegree, maxDegree));
67351a454edSToby Isaac   }
67451a454edSToby Isaac   PetscFunctionReturn(0);
67551a454edSToby Isaac }
67651a454edSToby Isaac 
677a2aba86cSMatthew G. Knepley PetscErrorCode DMFieldGetFVQuadrature_Internal(DMField field, IS pointIS, PetscQuadrature *quad)
678a2aba86cSMatthew G. Knepley {
679a2aba86cSMatthew G. Knepley   DM              dm = field->dm;
680a2aba86cSMatthew G. Knepley   const PetscInt *points;
681a2aba86cSMatthew G. Knepley   DMPolytopeType  ct;
682a2aba86cSMatthew G. Knepley   PetscInt        dim, n;
683a2aba86cSMatthew G. Knepley   PetscBool       isplex;
684a2aba86cSMatthew G. Knepley 
685a2aba86cSMatthew G. Knepley   PetscFunctionBegin;
686*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isplex));
687*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(pointIS, &n));
688a2aba86cSMatthew G. Knepley   if (isplex && n) {
689*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDimension(dm, &dim));
690*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(pointIS, &points));
691*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCellType(dm, points[0], &ct));
692a2aba86cSMatthew G. Knepley     switch (ct) {
693a2aba86cSMatthew G. Knepley       case DM_POLYTOPE_TRIANGLE:
694a2aba86cSMatthew G. Knepley       case DM_POLYTOPE_TETRAHEDRON:
695*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, 1, -1.0, 1.0, quad));break;
696*5f80ce2aSJacob Faibussowitsch       default: CHKERRQ(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
697a2aba86cSMatthew G. Knepley     }
698*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(pointIS, &points));
699a2aba86cSMatthew G. Knepley   } else {
700*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldCreateDefaultQuadrature(field, pointIS, quad));
701a2aba86cSMatthew G. Knepley   }
702a2aba86cSMatthew G. Knepley   PetscFunctionReturn(0);
703a2aba86cSMatthew G. Knepley }
704a2aba86cSMatthew G. Knepley 
70551a454edSToby Isaac static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad)
70651a454edSToby Isaac {
7071fdf1f07SMatthew G. Knepley   PetscInt       h, dim, imax, imin, cellHeight;
70851a454edSToby Isaac   DM             dm;
70951a454edSToby Isaac   DMField_DS     *dsfield;
71051a454edSToby Isaac   PetscObject    disc;
71151a454edSToby Isaac   PetscFE        fe;
71251a454edSToby Isaac   PetscClassId   id;
71351a454edSToby Isaac 
71451a454edSToby Isaac   PetscFunctionBegin;
71551a454edSToby Isaac   dm = field->dm;
71651a454edSToby Isaac   dsfield = (DMField_DS *) field->data;
717*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetMinMax(pointIS,&imin,&imax));
718*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
71951a454edSToby Isaac   for (h = 0; h <= dim; h++) {
72051a454edSToby Isaac     PetscInt hStart, hEnd;
72151a454edSToby Isaac 
722*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm,h,&hStart,&hEnd));
7236c041b70SSatish Balay     if (imax >= hStart && imin < hEnd) break;
72451a454edSToby Isaac   }
725*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetVTKCellHeight(dm, &cellHeight));
7261fdf1f07SMatthew G. Knepley   h -= cellHeight;
72751a454edSToby Isaac   *quad = NULL;
728f99c8401SToby Isaac   if (h < dsfield->height) {
729*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldDSGetHeightDisc(field,h,&disc));
730*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(disc,&id));
73151a454edSToby Isaac     if (id != PETSCFE_CLASSID) PetscFunctionReturn(0);
73251a454edSToby Isaac     fe = (PetscFE) disc;
733*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetQuadrature(fe,quad));
734*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)*quad));
73551a454edSToby Isaac   }
73651a454edSToby Isaac   PetscFunctionReturn(0);
73751a454edSToby Isaac }
73851a454edSToby Isaac 
7390145028aSToby Isaac static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom)
7400145028aSToby Isaac {
7410145028aSToby Isaac   const PetscInt *points;
7420145028aSToby Isaac   PetscInt        p, dim, dE, numFaces, Nq;
743b7260050SToby Isaac   PetscInt        maxDegree;
7440145028aSToby Isaac   DMLabel         depthLabel;
7450145028aSToby Isaac   IS              cellIS;
7460145028aSToby Isaac   DM              dm = field->dm;
7470145028aSToby Isaac 
7480145028aSToby Isaac   PetscFunctionBegin;
7490145028aSToby Isaac   dim = geom->dim;
7500145028aSToby Isaac   dE  = geom->dimEmbed;
751*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthLabel(dm, &depthLabel));
752*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS));
753*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldGetDegree(field,cellIS,NULL,&maxDegree));
754*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(pointIS, &points));
7550145028aSToby Isaac   numFaces = geom->numCells;
7560145028aSToby Isaac   Nq = geom->numPoints;
757f15274beSMatthew Knepley   /* First, set local faces and flip normals so that they are outward for the first supporting cell */
758f15274beSMatthew Knepley   for (p = 0; p < numFaces; p++) {
759f15274beSMatthew Knepley     PetscInt        point = points[p];
760f15274beSMatthew Knepley     PetscInt        suppSize, s, coneSize, c, numChildren;
761f15274beSMatthew Knepley     const PetscInt *supp, *cone, *ornt;
762f15274beSMatthew Knepley 
763*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
7642c71b3e2SJacob Faibussowitsch     PetscCheckFalse(numChildren,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
765*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetSupportSize(dm, point, &suppSize));
7662c71b3e2SJacob Faibussowitsch     PetscCheckFalse(suppSize > 2,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2", point, suppSize);
767f15274beSMatthew Knepley     if (!suppSize) continue;
768*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetSupport(dm, point, &supp));
769f15274beSMatthew Knepley     for (s = 0; s < suppSize; ++s) {
770*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetConeSize(dm, supp[s], &coneSize));
771*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetCone(dm, supp[s], &cone));
772f15274beSMatthew Knepley       for (c = 0; c < coneSize; ++c) if (cone[c] == point) break;
7732c71b3e2SJacob Faibussowitsch       PetscCheckFalse(c == coneSize,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid connectivity: point %D not found in cone of support point %D", point, supp[s]);
774f15274beSMatthew Knepley       geom->face[p][s] = c;
775f15274beSMatthew Knepley     }
776*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeOrientation(dm, supp[0], &ornt));
777f15274beSMatthew Knepley     if (ornt[geom->face[p][0]] < 0) {
778f15274beSMatthew Knepley       PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d;
779f15274beSMatthew Knepley 
780f15274beSMatthew Knepley       for (q = 0; q < Np; ++q) for (d = 0; d < dE; ++d) geom->n[(p*Np + q)*dE + d] = -geom->n[(p*Np + q)*dE + d];
781f15274beSMatthew Knepley     }
782f15274beSMatthew Knepley   }
783b7260050SToby Isaac   if (maxDegree <= 1) {
7840145028aSToby Isaac     PetscInt        numCells, offset, *cells;
7850145028aSToby Isaac     PetscFEGeom     *cellGeom;
7860145028aSToby Isaac     IS              suppIS;
7870145028aSToby Isaac 
7880145028aSToby Isaac     for (p = 0, numCells = 0; p < numFaces; p++) {
7890145028aSToby Isaac       PetscInt        point = points[p];
7900145028aSToby Isaac       PetscInt        numSupp, numChildren;
7910145028aSToby Isaac 
792*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
7932c71b3e2SJacob Faibussowitsch       PetscCheckFalse(numChildren,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
794*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm, point,&numSupp));
7952c71b3e2SJacob Faibussowitsch       PetscCheckFalse(numSupp > 2,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2", point, numSupp);
7960145028aSToby Isaac       numCells += numSupp;
7970145028aSToby Isaac     }
798*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(numCells, &cells));
7990145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8000145028aSToby Isaac       PetscInt        point = points[p];
8010145028aSToby Isaac       PetscInt        numSupp, s;
8020145028aSToby Isaac       const PetscInt *supp;
8030145028aSToby Isaac 
804*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm, point,&numSupp));
805*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupport(dm, point, &supp));
8060145028aSToby Isaac       for (s = 0; s < numSupp; s++, offset++) {
8070145028aSToby Isaac         cells[offset] = supp[s];
8080145028aSToby Isaac       }
8090145028aSToby Isaac     }
810*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS));
811*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldCreateFEGeom(field,suppIS,quad,PETSC_FALSE,&cellGeom));
8120145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8130145028aSToby Isaac       PetscInt        point = points[p];
8140145028aSToby Isaac       PetscInt        numSupp, s, q;
8150145028aSToby Isaac       const PetscInt *supp;
8160145028aSToby Isaac 
817*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm, point,&numSupp));
818*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupport(dm, point, &supp));
8190145028aSToby Isaac       for (s = 0; s < numSupp; s++, offset++) {
8200145028aSToby Isaac         for (q = 0; q < Nq * dE * dE; q++) {
82141418987SMatthew G. Knepley           geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
8220145028aSToby Isaac           geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
8230145028aSToby Isaac         }
82441418987SMatthew G. Knepley         for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
8250145028aSToby Isaac       }
8260145028aSToby Isaac     }
827*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGeomDestroy(&cellGeom));
828*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&suppIS));
829*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(cells));
8300145028aSToby Isaac   } else {
8310145028aSToby Isaac     PetscObject          faceDisc, cellDisc;
8320145028aSToby Isaac     PetscClassId         faceId, cellId;
8330145028aSToby Isaac     PetscDualSpace       dsp;
8340145028aSToby Isaac     DM                   K;
835ba2698f1SMatthew G. Knepley     DMPolytopeType       ct;
8360145028aSToby Isaac     PetscInt           (*co)[2][3];
8370145028aSToby Isaac     PetscInt             coneSize;
8380145028aSToby Isaac     PetscInt           **counts;
8390145028aSToby Isaac     PetscInt             f, i, o, q, s;
8400145028aSToby Isaac     const PetscInt      *coneK;
841ba2698f1SMatthew G. Knepley     PetscInt             eStart, minOrient, maxOrient, numOrient;
8420145028aSToby Isaac     PetscInt            *orients;
8430145028aSToby Isaac     PetscReal          **orientPoints;
8440145028aSToby Isaac     PetscReal           *cellPoints;
8450145028aSToby Isaac     PetscReal           *dummyWeights;
8460145028aSToby Isaac     PetscQuadrature      cellQuad = NULL;
8470145028aSToby Isaac 
848*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldDSGetHeightDisc(field, 1, &faceDisc));
849*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMFieldDSGetHeightDisc(field, 0, &cellDisc));
850*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(faceDisc,&faceId));
851*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(cellDisc,&cellId));
8522c71b3e2SJacob Faibussowitsch     PetscCheckFalse(faceId != PETSCFE_CLASSID || cellId != PETSCFE_CLASSID,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported");
853*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetDualSpace((PetscFE)cellDisc, &dsp));
854*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDualSpaceGetDM(dsp, &K));
855*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(K, 1, &eStart, NULL));
856*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCellType(K, eStart, &ct));
857*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeSize(K,0,&coneSize));
858*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCone(K,0,&coneK));
859*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(numFaces, &co, coneSize, &counts));
860*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(dE*Nq, &cellPoints));
861*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(Nq, &dummyWeights));
862*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad));
863*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights));
8640145028aSToby Isaac     minOrient = PETSC_MAX_INT;
8650145028aSToby Isaac     maxOrient = PETSC_MIN_INT;
8660145028aSToby Isaac     for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */
8670145028aSToby Isaac       PetscInt        point = points[p];
8680145028aSToby Isaac       PetscInt        numSupp, numChildren;
8690145028aSToby Isaac       const PetscInt *supp;
8700145028aSToby Isaac 
871*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
8722c71b3e2SJacob Faibussowitsch       PetscCheckFalse(numChildren,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
873*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm, point,&numSupp));
8742c71b3e2SJacob Faibussowitsch       PetscCheckFalse(numSupp > 2,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2", point, numSupp);
875*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupport(dm, point, &supp));
8760145028aSToby Isaac       for (s = 0; s < numSupp; s++) {
8770145028aSToby Isaac         PetscInt        cell = supp[s];
8780145028aSToby Isaac         PetscInt        numCone;
8790145028aSToby Isaac         const PetscInt *cone, *orient;
8800145028aSToby Isaac 
881*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetConeSize(dm, cell, &numCone));
8822c71b3e2SJacob Faibussowitsch         PetscCheckFalse(numCone != coneSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element");
883*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetCone(dm, cell, &cone));
884*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetConeOrientation(dm, cell, &orient));
8850145028aSToby Isaac         for (f = 0; f < coneSize; f++) {
8860145028aSToby Isaac           if (cone[f] == point) break;
8870145028aSToby Isaac         }
8880145028aSToby Isaac         co[p][s][0] = f;
8890145028aSToby Isaac         co[p][s][1] = orient[f];
8900145028aSToby Isaac         co[p][s][2] = cell;
8910145028aSToby Isaac         minOrient = PetscMin(minOrient, orient[f]);
892cd93b0e1SMatthew G. Knepley         maxOrient = PetscMax(maxOrient, orient[f]);
8930145028aSToby Isaac       }
8940145028aSToby Isaac       for (; s < 2; s++) {
8950145028aSToby Isaac         co[p][s][0] = -1;
8960145028aSToby Isaac         co[p][s][1] = -1;
8970145028aSToby Isaac         co[p][s][2] = -1;
8980145028aSToby Isaac       }
8990145028aSToby Isaac     }
9000145028aSToby Isaac     numOrient = maxOrient + 1 - minOrient;
901*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCone(K,0,&coneK));
9020145028aSToby Isaac     /* count all (face,orientation) doubles that appear */
903*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc2(numOrient,&orients,numOrient,&orientPoints));
904*5f80ce2aSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) CHKERRQ(PetscCalloc1(numOrient+1, &counts[f]));
9050145028aSToby Isaac     for (p = 0; p < numFaces; p++) {
9060145028aSToby Isaac       for (s = 0; s < 2; s++) {
9070145028aSToby Isaac         if (co[p][s][0] >= 0) {
9080145028aSToby Isaac           counts[co[p][s][0]][co[p][s][1] - minOrient]++;
9090145028aSToby Isaac           orients[co[p][s][1] - minOrient]++;
9100145028aSToby Isaac         }
9110145028aSToby Isaac       }
9120145028aSToby Isaac     }
9130145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
9140145028aSToby Isaac       if (orients[o]) {
9150145028aSToby Isaac         PetscInt orient = o + minOrient;
9160145028aSToby Isaac         PetscInt q;
9170145028aSToby Isaac 
918*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(Nq * dim, &orientPoints[o]));
9190145028aSToby Isaac         /* rotate the quadrature points appropriately */
920ba2698f1SMatthew G. Knepley         switch (ct) {
921ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_POINT: break;
922ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_SEGMENT:
9230145028aSToby Isaac           if (orient == -2 || orient == 1) {
9240145028aSToby Isaac             for (q = 0; q < Nq; q++) {
9250145028aSToby Isaac               orientPoints[o][q] = -geom->xi[q];
9260145028aSToby Isaac             }
9270145028aSToby Isaac           } else {
9280145028aSToby Isaac             for (q = 0; q < Nq; q++) {
9290145028aSToby Isaac               orientPoints[o][q] = geom->xi[q];
9300145028aSToby Isaac             }
9310145028aSToby Isaac           }
9320145028aSToby Isaac           break;
933ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:
9340145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9350145028aSToby Isaac             PetscReal lambda[3];
9360145028aSToby Isaac             PetscReal lambdao[3];
9370145028aSToby Isaac 
9380145028aSToby Isaac             /* convert to barycentric */
9390145028aSToby Isaac             lambda[0] = - (geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.;
9400145028aSToby Isaac             lambda[1] = (geom->xi[2 * q] + 1.) / 2.;
9410145028aSToby Isaac             lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.;
9420145028aSToby Isaac             if (orient >= 0) {
9430145028aSToby Isaac               for (i = 0; i < 3; i++) {
9440145028aSToby Isaac                 lambdao[i] = lambda[(orient + i) % 3];
9450145028aSToby Isaac               }
9460145028aSToby Isaac             } else {
9470145028aSToby Isaac               for (i = 0; i < 3; i++) {
9480145028aSToby Isaac                 lambdao[i] = lambda[(-(orient + i) + 3) % 3];
9490145028aSToby Isaac               }
9500145028aSToby Isaac             }
9510145028aSToby Isaac             /* convert to coordinates */
9520145028aSToby Isaac             orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1];
9530145028aSToby Isaac             orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2];
9540145028aSToby Isaac           }
9550145028aSToby Isaac           break;
956ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL:
9570145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9580145028aSToby Isaac             PetscReal xi[2], xio[2];
9590145028aSToby Isaac             PetscInt oabs = (orient >= 0) ? orient : -(orient + 1);
9600145028aSToby Isaac 
9610145028aSToby Isaac             xi[0] = geom->xi[2 * q];
9620145028aSToby Isaac             xi[1] = geom->xi[2 * q + 1];
9630145028aSToby Isaac             switch (oabs) {
9640145028aSToby Isaac             case 1:
9650145028aSToby Isaac               xio[0] = xi[1];
9660145028aSToby Isaac               xio[1] = -xi[0];
9670145028aSToby Isaac               break;
9680145028aSToby Isaac             case 2:
9690145028aSToby Isaac               xio[0] = -xi[0];
9700145028aSToby Isaac               xio[1] = -xi[1];
9710145028aSToby Isaac             case 3:
9720145028aSToby Isaac               xio[0] = -xi[1];
9730145028aSToby Isaac               xio[1] = xi[0];
974aa1371cbSToby Isaac             case 0:
975aa1371cbSToby Isaac             default:
976aa1371cbSToby Isaac               xio[0] = xi[0];
977aa1371cbSToby Isaac               xio[1] = xi[1];
978aa1371cbSToby Isaac               break;
9790145028aSToby Isaac             }
9800145028aSToby Isaac             if (orient < 0) {
9810145028aSToby Isaac               xio[0] = -xio[0];
9820145028aSToby Isaac             }
9830145028aSToby Isaac             orientPoints[o][2 * q + 0] = xio[0];
9840145028aSToby Isaac             orientPoints[o][2 * q + 1] = xio[1];
9850145028aSToby Isaac           }
9860145028aSToby Isaac           break;
98798921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell type %s not yet supported", DMPolytopeTypes[ct]);
9880145028aSToby Isaac         }
9890145028aSToby Isaac       }
9900145028aSToby Isaac     }
9910145028aSToby Isaac     for (f = 0; f < coneSize; f++) {
9920145028aSToby Isaac       PetscInt face = coneK[f];
9934d1a973fSToby Isaac       PetscReal v0[3];
9942d89661fSMatthew G. Knepley       PetscReal J[9], detJ;
9950145028aSToby Isaac       PetscInt numCells, offset;
9960145028aSToby Isaac       PetscInt *cells;
9970145028aSToby Isaac       IS suppIS;
9980145028aSToby Isaac 
999*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ));
10000145028aSToby Isaac       for (o = 0; o <= numOrient; o++) {
10010145028aSToby Isaac         PetscFEGeom *cellGeom;
10020145028aSToby Isaac 
10030145028aSToby Isaac         if (!counts[f][o]) continue;
10040145028aSToby Isaac         /* If this (face,orientation) double appears,
10050145028aSToby Isaac          * convert the face quadrature points into volume quadrature points */
10060145028aSToby Isaac         for (q = 0; q < Nq; q++) {
10070145028aSToby Isaac           PetscReal xi0[3] = {-1., -1., -1.};
10080145028aSToby Isaac 
10090145028aSToby Isaac           CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]);
10100145028aSToby Isaac         }
10110145028aSToby Isaac         for (p = 0, numCells = 0; p < numFaces; p++) {
10120145028aSToby Isaac           for (s = 0; s < 2; s++) {
10130145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++;
10140145028aSToby Isaac           }
10150145028aSToby Isaac         }
1016*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(numCells, &cells));
10170145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10180145028aSToby Isaac           for (s = 0; s < 2; s++) {
10190145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
10200145028aSToby Isaac               cells[offset++] = co[p][s][2];
10210145028aSToby Isaac             }
10220145028aSToby Isaac           }
10230145028aSToby Isaac         }
1024*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS));
1025*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom));
10260145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10270145028aSToby Isaac           for (s = 0; s < 2; s++) {
10280145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
10290145028aSToby Isaac               for (q = 0; q < Nq * dE * dE; q++) {
103041418987SMatthew G. Knepley                 geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
10310145028aSToby Isaac                 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
10320145028aSToby Isaac               }
103341418987SMatthew G. Knepley               for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
10340145028aSToby Isaac               offset++;
10350145028aSToby Isaac             }
10360145028aSToby Isaac           }
10370145028aSToby Isaac         }
1038*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFEGeomDestroy(&cellGeom));
1039*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&suppIS));
1040*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(cells));
10410145028aSToby Isaac       }
10420145028aSToby Isaac     }
10430145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
10440145028aSToby Isaac       if (orients[o]) {
1045*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(orientPoints[o]));
10460145028aSToby Isaac       }
10470145028aSToby Isaac     }
1048*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(orients,orientPoints));
1049*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscQuadratureDestroy(&cellQuad));
1050*5f80ce2aSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) CHKERRQ(PetscFree(counts[f]));
1051*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(co,counts));
10520145028aSToby Isaac   }
1053*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(pointIS, &points));
1054*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&cellIS));
10550145028aSToby Isaac   PetscFunctionReturn(0);
10560145028aSToby Isaac }
10570145028aSToby Isaac 
105851a454edSToby Isaac static PetscErrorCode DMFieldInitialize_DS(DMField field)
105951a454edSToby Isaac {
106051a454edSToby Isaac   PetscFunctionBegin;
106151a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DS;
106251a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DS;
106351a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DS;
1064805e7170SToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DS;
1065b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DS;
106651a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS;
106751a454edSToby Isaac   field->ops->view                    = DMFieldView_DS;
10680145028aSToby Isaac   field->ops->computeFaceData         = DMFieldComputeFaceData_DS;
106951a454edSToby Isaac   PetscFunctionReturn(0);
107051a454edSToby Isaac }
107151a454edSToby Isaac 
107251a454edSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field)
107351a454edSToby Isaac {
107451a454edSToby Isaac   DMField_DS     *dsfield;
107551a454edSToby Isaac 
107651a454edSToby Isaac   PetscFunctionBegin;
1077*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(field,&dsfield));
107851a454edSToby Isaac   field->data = dsfield;
1079*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldInitialize_DS(field));
108051a454edSToby Isaac   PetscFunctionReturn(0);
108151a454edSToby Isaac }
108251a454edSToby Isaac 
108351a454edSToby Isaac PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec,DMField *field)
108451a454edSToby Isaac {
108551a454edSToby Isaac   DMField        b;
108651a454edSToby Isaac   DMField_DS     *dsfield;
108751a454edSToby Isaac   PetscObject    disc = NULL;
108851a454edSToby Isaac   PetscBool      isContainer = PETSC_FALSE;
108951a454edSToby Isaac   PetscClassId   id = -1;
109051a454edSToby Isaac   PetscInt       numComponents = -1, dsNumFields;
109151a454edSToby Isaac   PetscSection   section;
109251a454edSToby Isaac 
109351a454edSToby Isaac   PetscFunctionBegin;
1094*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm,&section));
1095*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetFieldComponents(section,fieldNum,&numComponents));
1096*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNumFields(dm,&dsNumFields));
1097*5f80ce2aSJacob Faibussowitsch   if (dsNumFields) CHKERRQ(DMGetField(dm,fieldNum,NULL,&disc));
109851a454edSToby Isaac   if (disc) {
1099*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(disc,&id));
110051a454edSToby Isaac     isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE;
110151a454edSToby Isaac   }
110251a454edSToby Isaac   if (!disc || isContainer) {
110351a454edSToby Isaac     MPI_Comm       comm = PetscObjectComm((PetscObject) dm);
110451a454edSToby Isaac     PetscFE        fe;
11052df84da0SMatthew G. Knepley     DMPolytopeType ct, locct = DM_POLYTOPE_UNKNOWN;
11062df84da0SMatthew G. Knepley     PetscInt       dim, cStart, cEnd, cellHeight;
110751a454edSToby Isaac 
1108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetVTKCellHeight(dm, &cellHeight));
1109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDimension(dm, &dim));
1110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
1111*5f80ce2aSJacob Faibussowitsch     if (cEnd > cStart) CHKERRQ(DMPlexGetCellType(dm, cStart, &locct));
1112*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Allreduce(&locct, &ct, 1, MPI_INT, MPI_MIN, comm));
1113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, numComponents, ct, 1, PETSC_DETERMINE, &fe));
1114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEViewFromOptions(fe, NULL, "-field_fe_view"));
111551a454edSToby Isaac     disc = (PetscObject) fe;
111651a454edSToby Isaac   } else {
1117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference(disc));
111851a454edSToby Isaac   }
1119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetClassId(disc,&id));
112051a454edSToby Isaac   if (id == PETSCFE_CLASSID) {
112151a454edSToby Isaac     PetscFE fe = (PetscFE) disc;
112251a454edSToby Isaac 
1123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEGetNumComponents(fe,&numComponents));
11242da392ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented");
1125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldCreate(dm,numComponents,DMFIELD_VERTEX,&b));
1126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMFieldSetType(b,DMFIELDDS));
112751a454edSToby Isaac   dsfield = (DMField_DS *) b->data;
112851a454edSToby Isaac   dsfield->fieldNum = fieldNum;
1129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dsfield->height));
1130f99c8401SToby Isaac   dsfield->height++;
1131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(dsfield->height,&dsfield->disc));
113251a454edSToby Isaac   dsfield->disc[0] = disc;
1133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)vec));
113451a454edSToby Isaac   dsfield->vec = vec;
113551a454edSToby Isaac   *field = b;
113651a454edSToby Isaac 
113751a454edSToby Isaac   PetscFunctionReturn(0);
113851a454edSToby Isaac }
1139