xref: /petsc/src/dm/field/impls/ds/dmfieldds.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
79371c9d4SSatish Balay typedef struct _n_DMField_DS {
851a454edSToby Isaac   PetscBool    multifieldVec;
96858538eSMatthew G. Knepley   PetscInt     height;   /* Point height at which we want values and number of discretizations */
106858538eSMatthew G. Knepley   PetscInt     fieldNum; /* Number in DS of field which we evaluate */
116858538eSMatthew G. Knepley   PetscObject *disc;     /* Discretizations of this field at each height */
126858538eSMatthew G. Knepley   Vec          vec;      /* Field values */
136858538eSMatthew G. Knepley   DM           dmDG;     /* DM for the DG values */
146858538eSMatthew G. Knepley   PetscObject *discDG;   /* DG Discretizations of this field at each height */
156858538eSMatthew G. Knepley   Vec          vecDG;    /* DG Field values */
166858538eSMatthew G. Knepley } DMField_DS;
1751a454edSToby Isaac 
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldDestroy_DS(DMField field)
19d71ae5a4SJacob Faibussowitsch {
206858538eSMatthew G. Knepley   DMField_DS *dsfield = (DMField_DS *)field->data;
2151a454edSToby Isaac   PetscInt    i;
2251a454edSToby Isaac 
2351a454edSToby Isaac   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dsfield->vec));
256858538eSMatthew G. Knepley   for (i = 0; i < dsfield->height; i++) PetscCall(PetscObjectDereference(dsfield->disc[i]));
269566063dSJacob Faibussowitsch   PetscCall(PetscFree(dsfield->disc));
276858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dsfield->vecDG));
289371c9d4SSatish Balay   if (dsfield->discDG)
299371c9d4SSatish Balay     for (i = 0; i < dsfield->height; i++) PetscCall(PetscObjectDereference(dsfield->discDG[i]));
306858538eSMatthew G. Knepley   PetscCall(PetscFree(dsfield->discDG));
319566063dSJacob Faibussowitsch   PetscCall(PetscFree(dsfield));
32*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3351a454edSToby Isaac }
3451a454edSToby Isaac 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldView_DS(DMField field, PetscViewer viewer)
36d71ae5a4SJacob Faibussowitsch {
3751a454edSToby Isaac   DMField_DS *dsfield = (DMField_DS *)field->data;
3851a454edSToby Isaac   PetscObject disc;
396858538eSMatthew G. Knepley   PetscBool   iascii;
4051a454edSToby Isaac 
4151a454edSToby Isaac   PetscFunctionBegin;
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4351a454edSToby Isaac   disc = dsfield->disc[0];
4451a454edSToby Isaac   if (iascii) {
4563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PetscDS field %" PetscInt_FMT "\n", dsfield->fieldNum));
469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
479566063dSJacob Faibussowitsch     PetscCall(PetscObjectView(disc, viewer));
489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
4951a454edSToby Isaac   }
509566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
5128b400f6SJacob Faibussowitsch   PetscCheck(!dsfield->multifieldVec, PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "View of subfield not implemented yet");
529566063dSJacob Faibussowitsch   PetscCall(VecView(dsfield->vec, viewer));
539566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
54*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5551a454edSToby Isaac }
5651a454edSToby Isaac 
57d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject discList[], PetscObject *disc)
58d71ae5a4SJacob Faibussowitsch {
5951a454edSToby Isaac   PetscFunctionBegin;
606858538eSMatthew G. Knepley   if (!discList[height]) {
6151a454edSToby Isaac     PetscClassId id;
6251a454edSToby Isaac 
636858538eSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(discList[0], &id));
646858538eSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateHeightTrace((PetscFE)discList[0], height, (PetscFE *)&discList[height]));
6551a454edSToby Isaac   }
666858538eSMatthew G. Knepley   *disc = discList[height];
67*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6851a454edSToby Isaac }
6951a454edSToby Isaac 
7062bd480fSMatthew G. Knepley /* y[m,c] = A[m,n,c] . b[n] */
7151a454edSToby Isaac #define DMFieldDSdot(y, A, b, m, n, c, cast) \
7251a454edSToby Isaac   do { \
7351a454edSToby Isaac     PetscInt _i, _j, _k; \
7451a454edSToby Isaac     for (_i = 0; _i < (m); _i++) { \
75ad540459SPierre Jolivet       for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] = 0.; \
7651a454edSToby Isaac       for (_j = 0; _j < (n); _j++) { \
77ad540459SPierre Jolivet         for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \
7851a454edSToby Isaac       } \
7951a454edSToby Isaac     } \
8051a454edSToby Isaac   } while (0)
8151a454edSToby Isaac 
826858538eSMatthew G. Knepley /*
836858538eSMatthew G. Knepley   Since this is used for coordinates, we need to allow for the possibility that values come from multiple sections/Vecs, so that we can have DG version of the coordinates for periodicity. This reproduces DMPlexGetCellCoordinates_Internal().
846858538eSMatthew G. Knepley */
85d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[])
86d71ae5a4SJacob Faibussowitsch {
876858538eSMatthew G. Knepley   DMField_DS        *dsfield = (DMField_DS *)field->data;
886858538eSMatthew G. Knepley   DM                 fdm     = dsfield->dmDG;
896858538eSMatthew G. Knepley   PetscSection       s       = NULL;
906858538eSMatthew G. Knepley   const PetscScalar *cvalues;
916858538eSMatthew G. Knepley   PetscInt           pStart, pEnd;
926858538eSMatthew G. Knepley 
936858538eSMatthew G. Knepley   PetscFunctionBeginHot;
946858538eSMatthew G. Knepley   *isDG   = PETSC_FALSE;
956858538eSMatthew G. Knepley   *Nc     = 0;
966858538eSMatthew G. Knepley   *array  = NULL;
976858538eSMatthew G. Knepley   *values = NULL;
986858538eSMatthew G. Knepley   /* Check for cellwise section */
996858538eSMatthew G. Knepley   if (fdm) PetscCall(DMGetLocalSection(fdm, &s));
1006858538eSMatthew G. Knepley   if (!s) goto cg;
1016858538eSMatthew G. Knepley   /* Check that the cell exists in the cellwise section */
1026858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1036858538eSMatthew G. Knepley   if (cell < pStart || cell >= pEnd) goto cg;
1046858538eSMatthew G. Knepley   /* Check for cellwise coordinates for this cell */
1056858538eSMatthew G. Knepley   PetscCall(PetscSectionGetDof(s, cell, Nc));
1066858538eSMatthew G. Knepley   if (!*Nc) goto cg;
1076858538eSMatthew G. Knepley   /* Check for cellwise coordinates */
1086858538eSMatthew G. Knepley   if (!dsfield->vecDG) goto cg;
1096858538eSMatthew G. Knepley   /* Get cellwise coordinates */
1106858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(dsfield->vecDG, array));
1116858538eSMatthew G. Knepley   PetscCall(DMPlexPointLocalRead(fdm, cell, *array, &cvalues));
1126858538eSMatthew G. Knepley   PetscCall(DMGetWorkArray(fdm, *Nc, MPIU_SCALAR, values));
1136858538eSMatthew G. Knepley   PetscCall(PetscArraycpy(*values, cvalues, *Nc));
1146858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(dsfield->vecDG, array));
1156858538eSMatthew G. Knepley   *isDG = PETSC_TRUE;
116*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1176858538eSMatthew G. Knepley cg:
1186858538eSMatthew G. Knepley   /* Use continuous values */
1196858538eSMatthew G. Knepley   PetscCall(DMFieldGetDM(field, &fdm));
1206858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(fdm, &s));
1216267c18cSMatthew G. Knepley   PetscCall(PetscSectionGetField(s, dsfield->fieldNum, &s));
1226858538eSMatthew G. Knepley   PetscCall(DMPlexVecGetClosure(fdm, s, dsfield->vec, cell, Nc, values));
123*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1246858538eSMatthew G. Knepley }
1256858538eSMatthew G. Knepley 
126d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldRestoreClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[])
127d71ae5a4SJacob Faibussowitsch {
1286858538eSMatthew G. Knepley   DMField_DS  *dsfield = (DMField_DS *)field->data;
1296858538eSMatthew G. Knepley   DM           fdm;
1306858538eSMatthew G. Knepley   PetscSection s;
1316858538eSMatthew G. Knepley 
1326858538eSMatthew G. Knepley   PetscFunctionBeginHot;
1336858538eSMatthew G. Knepley   if (*isDG) {
1346858538eSMatthew G. Knepley     PetscCall(DMRestoreWorkArray(dsfield->dmDG, *Nc, MPIU_SCALAR, values));
1356858538eSMatthew G. Knepley   } else {
1366858538eSMatthew G. Knepley     PetscCall(DMFieldGetDM(field, &fdm));
1376858538eSMatthew G. Knepley     PetscCall(DMGetLocalSection(fdm, &s));
1386858538eSMatthew G. Knepley     PetscCall(DMPlexVecRestoreClosure(fdm, s, dsfield->vec, cell, Nc, (PetscScalar **)values));
1396858538eSMatthew G. Knepley   }
140*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1416858538eSMatthew G. Knepley }
1426858538eSMatthew G. Knepley 
143ef0bb6c7SMatthew G. Knepley /* TODO: Reorganize interface so that I can reuse a tabulation rather than mallocing each time */
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H)
145d71ae5a4SJacob Faibussowitsch {
14651a454edSToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
14751a454edSToby Isaac   DM               dm;
14851a454edSToby Isaac   PetscObject      disc;
14951a454edSToby Isaac   PetscClassId     classid;
15051a454edSToby Isaac   PetscInt         nq, nc, dim, meshDim, numCells;
15151a454edSToby Isaac   PetscSection     section;
15251a454edSToby Isaac   const PetscReal *qpoints;
15351a454edSToby Isaac   PetscBool        isStride;
15451a454edSToby Isaac   const PetscInt  *points = NULL;
15551a454edSToby Isaac   PetscInt         sfirst = -1, stride = -1;
15651a454edSToby Isaac 
15751a454edSToby Isaac   PetscFunctionBeginHot;
15851a454edSToby Isaac   dm = field->dm;
15951a454edSToby Isaac   nc = field->numComponents;
1609566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &nq, &qpoints, NULL));
1616858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, dsfield->height - 1 - dim, dsfield->disc, &disc));
1629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &meshDim));
1639566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1649566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetField(section, dsfield->fieldNum, &section));
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &classid));
16651a454edSToby Isaac   /* TODO: batch */
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
1689566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
1691baa6e33SBarry Smith   if (isStride) PetscCall(ISStrideGetInfo(pointIS, &sfirst, &stride));
1701baa6e33SBarry Smith   else PetscCall(ISGetIndices(pointIS, &points));
17151a454edSToby Isaac   if (classid == PETSCFE_CLASSID) {
17251a454edSToby Isaac     PetscFE         fe = (PetscFE)disc;
17351a454edSToby Isaac     PetscInt        feDim, i;
174ef0bb6c7SMatthew G. Knepley     PetscInt        K = H ? 2 : (D ? 1 : (B ? 0 : -1));
175ef0bb6c7SMatthew G. Knepley     PetscTabulation T;
17651a454edSToby Isaac 
1779566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDimension(fe, &feDim));
1789566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fe, 1, nq, qpoints, K, &T));
17951a454edSToby Isaac     for (i = 0; i < numCells; i++) {
18051a454edSToby Isaac       PetscInt           c = isStride ? (sfirst + i * stride) : points[i];
1812710589cSToby Isaac       PetscInt           closureSize;
1826858538eSMatthew G. Knepley       const PetscScalar *array;
1832710589cSToby Isaac       PetscScalar       *elem = NULL;
1846858538eSMatthew G. Knepley       PetscBool          isDG;
18551a454edSToby Isaac 
1866858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
18751a454edSToby Isaac       if (B) {
18862bd480fSMatthew G. Knepley         /* field[c] = T[q,b,c] . coef[b], so v[c] = T[q,b,c] . coords[b] */
18951a454edSToby Isaac         if (type == PETSC_SCALAR) {
19051a454edSToby Isaac           PetscScalar *cB = &((PetscScalar *)B)[nc * nq * i];
19151a454edSToby Isaac 
192ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
19351a454edSToby Isaac         } else {
19451a454edSToby Isaac           PetscReal *cB = &((PetscReal *)B)[nc * nq * i];
19551a454edSToby Isaac 
196ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
19751a454edSToby Isaac         }
19851a454edSToby Isaac       }
19951a454edSToby Isaac       if (D) {
20051a454edSToby Isaac         if (type == PETSC_SCALAR) {
20151a454edSToby Isaac           PetscScalar *cD = &((PetscScalar *)D)[nc * nq * dim * i];
20251a454edSToby Isaac 
203ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
20451a454edSToby Isaac         } else {
20551a454edSToby Isaac           PetscReal *cD = &((PetscReal *)D)[nc * nq * dim * i];
20651a454edSToby Isaac 
207ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
20851a454edSToby Isaac         }
20951a454edSToby Isaac       }
21051a454edSToby Isaac       if (H) {
21151a454edSToby Isaac         if (type == PETSC_SCALAR) {
21251a454edSToby Isaac           PetscScalar *cH = &((PetscScalar *)H)[nc * nq * dim * dim * i];
21351a454edSToby Isaac 
214ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
21551a454edSToby Isaac         } else {
21651a454edSToby Isaac           PetscReal *cH = &((PetscReal *)H)[nc * nq * dim * dim * i];
21751a454edSToby Isaac 
218ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
21951a454edSToby Isaac         }
22051a454edSToby Isaac       }
2216858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
22251a454edSToby Isaac     }
2239566063dSJacob Faibussowitsch     PetscCall(PetscTabulationDestroy(&T));
2242da392ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented");
22548a46eb9SPierre Jolivet   if (!isStride) PetscCall(ISRestoreIndices(pointIS, &points));
226*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22751a454edSToby Isaac }
22851a454edSToby Isaac 
229d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
230d71ae5a4SJacob Faibussowitsch {
231a6216909SToby Isaac   DMField_DS        *dsfield = (DMField_DS *)field->data;
232a6216909SToby Isaac   PetscSF            cellSF  = NULL;
233a6216909SToby Isaac   const PetscSFNode *cells;
234a6216909SToby Isaac   PetscInt           c, nFound, numCells, feDim, nc;
235a6216909SToby Isaac   const PetscInt    *cellDegrees;
236a6216909SToby Isaac   const PetscScalar *pointsArray;
237a6216909SToby Isaac   PetscScalar       *cellPoints;
238a6216909SToby Isaac   PetscInt           gatherSize, gatherMax;
239a6216909SToby Isaac   PetscInt           dim, dimR, offset;
240a6216909SToby Isaac   MPI_Datatype       pointType;
241a6216909SToby Isaac   PetscObject        cellDisc;
242a6216909SToby Isaac   PetscFE            cellFE;
243a6216909SToby Isaac   PetscClassId       discID;
244a6216909SToby Isaac   PetscReal         *coordsReal, *coordsRef;
245a6216909SToby Isaac   PetscSection       section;
246a6216909SToby Isaac   PetscScalar       *cellBs = NULL, *cellDs = NULL, *cellHs = NULL;
247a6216909SToby Isaac   PetscReal         *cellBr = NULL, *cellDr = NULL, *cellHr = NULL;
2484cbe5319SToby Isaac   PetscReal         *v, *J, *invJ, *detJ;
249a6216909SToby Isaac 
250a6216909SToby Isaac   PetscFunctionBegin;
251a6216909SToby Isaac   nc = field->numComponents;
2529566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(field->dm, &section));
2536858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
2549566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(cellDisc, &discID));
25508401ef6SPierre Jolivet   PetscCheck(discID == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization type not supported");
256a6216909SToby Isaac   cellFE = (PetscFE)cellDisc;
2579566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDimension(cellFE, &feDim));
2589566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dim));
2599566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dimR));
2609566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF));
2619566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells));
262ad540459SPierre Jolivet   for (c = 0; c < nFound; c++) PetscCheck(cells[c].index >= 0, PetscObjectComm((PetscObject)points), PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " could not be located", c);
2639566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(cellSF, &cellDegrees));
2649566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(cellSF, &cellDegrees));
265a6216909SToby Isaac   for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) {
266a6216909SToby Isaac     gatherMax = PetscMax(gatherMax, cellDegrees[c]);
267a6216909SToby Isaac     gatherSize += cellDegrees[c];
268a6216909SToby Isaac   }
2699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(gatherSize * dim, &cellPoints, gatherMax * dim, &coordsReal, gatherMax * dimR, &coordsRef));
2709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(gatherMax * dimR, &v, gatherMax * dimR * dimR, &J, gatherMax * dimR * dimR, &invJ, gatherMax, &detJ));
2711baa6e33SBarry Smith   if (datatype == PETSC_SCALAR) PetscCall(PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs));
2721baa6e33SBarry Smith   else PetscCall(PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr));
273a6216909SToby Isaac 
2749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(dim, MPIU_SCALAR, &pointType));
2759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&pointType));
2769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points, &pointsArray));
2779566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints));
2789566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints));
2799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points, &pointsArray));
280a6216909SToby Isaac   for (c = 0, offset = 0; c < numCells; c++) {
281a6216909SToby Isaac     PetscInt nq = cellDegrees[c], p;
282a6216909SToby Isaac 
283a6216909SToby Isaac     if (nq) {
284ef0bb6c7SMatthew G. Knepley       PetscInt           K = H ? 2 : (D ? 1 : (B ? 0 : -1));
285ef0bb6c7SMatthew G. Knepley       PetscTabulation    T;
2866858538eSMatthew G. Knepley       PetscQuadrature    quad;
2876858538eSMatthew G. Knepley       const PetscScalar *array;
288a6216909SToby Isaac       PetscScalar       *elem = NULL;
2894cbe5319SToby Isaac       PetscReal         *quadPoints;
2906858538eSMatthew G. Knepley       PetscBool          isDG;
2916858538eSMatthew G. Knepley       PetscInt           closureSize, d, e, f, g;
292a6216909SToby Isaac 
2934d1a973fSToby Isaac       for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]);
2949566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef));
2959566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(cellFE, 1, nq, coordsRef, K, &T));
2969566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad));
2979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dimR * nq, &quadPoints));
2984cbe5319SToby Isaac       for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p];
2999566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL));
3009566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ));
3019566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
3026858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
303a6216909SToby Isaac       if (B) {
304a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
305a6216909SToby Isaac           PetscScalar *cB = &cellBs[nc * offset];
306a6216909SToby Isaac 
307ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
308a6216909SToby Isaac         } else {
309a6216909SToby Isaac           PetscReal *cB = &cellBr[nc * offset];
310a6216909SToby Isaac 
311ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
312a6216909SToby Isaac         }
313a6216909SToby Isaac       }
314a6216909SToby Isaac       if (D) {
315a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
316a6216909SToby Isaac           PetscScalar *cD = &cellDs[nc * dim * offset];
317a6216909SToby Isaac 
318ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
3194cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3204cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3214d1a973fSToby Isaac               PetscScalar vs[3];
3224d1a973fSToby Isaac 
3234cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3244d1a973fSToby Isaac                 vs[d] = 0.;
325ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3264cbe5319SToby Isaac               }
327ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = vs[d];
3284cbe5319SToby Isaac             }
3294cbe5319SToby Isaac           }
330a6216909SToby Isaac         } else {
331a6216909SToby Isaac           PetscReal *cD = &cellDr[nc * dim * offset];
332a6216909SToby Isaac 
333ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
3344cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3354cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3364cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3374cbe5319SToby Isaac                 v[d] = 0.;
338ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3394cbe5319SToby Isaac               }
340ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = v[d];
3414cbe5319SToby Isaac             }
3424cbe5319SToby Isaac           }
343a6216909SToby Isaac         }
344a6216909SToby Isaac       }
345a6216909SToby Isaac       if (H) {
346a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
347a6216909SToby Isaac           PetscScalar *cH = &cellHs[nc * dim * dim * offset];
348a6216909SToby Isaac 
349ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
3504cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3514cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3524d1a973fSToby Isaac               PetscScalar vs[3];
3534d1a973fSToby Isaac 
3544cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3554d1a973fSToby Isaac                 vs[d] = 0.;
356ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3574cbe5319SToby Isaac               }
358ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = vs[d];
3594cbe5319SToby Isaac             }
3604cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3614cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3624d1a973fSToby Isaac                 PetscScalar vs[3];
3634d1a973fSToby Isaac 
3644cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3654d1a973fSToby Isaac                   vs[d] = 0.;
366ad540459SPierre Jolivet                   for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3674cbe5319SToby Isaac                 }
368ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d];
3694cbe5319SToby Isaac               }
3704cbe5319SToby Isaac             }
3714cbe5319SToby Isaac           }
372a6216909SToby Isaac         } else {
373a6216909SToby Isaac           PetscReal *cH = &cellHr[nc * dim * dim * offset];
374a6216909SToby Isaac 
375ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
3764cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3774cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3784cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3794cbe5319SToby Isaac                 v[d] = 0.;
380ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3814cbe5319SToby Isaac               }
382ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = v[d];
3834cbe5319SToby Isaac             }
3844cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3854cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3864cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3874cbe5319SToby Isaac                   v[d] = 0.;
388ad540459SPierre Jolivet                   for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3894cbe5319SToby Isaac                 }
390ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = v[d];
3914cbe5319SToby Isaac               }
3924cbe5319SToby Isaac             }
3934cbe5319SToby Isaac           }
394a6216909SToby Isaac         }
395a6216909SToby Isaac       }
3966858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
3979566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&T));
398a6216909SToby Isaac     }
399a6216909SToby Isaac     offset += nq;
400a6216909SToby Isaac   }
401a6216909SToby Isaac   {
402a6216909SToby Isaac     MPI_Datatype origtype;
4034cbe5319SToby Isaac 
404a6216909SToby Isaac     if (datatype == PETSC_SCALAR) {
405a6216909SToby Isaac       origtype = MPIU_SCALAR;
406a6216909SToby Isaac     } else {
407a6216909SToby Isaac       origtype = MPIU_REAL;
408a6216909SToby Isaac     }
409a6216909SToby Isaac     if (B) {
410a6216909SToby Isaac       MPI_Datatype Btype;
411a6216909SToby Isaac 
4129566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc, origtype, &Btype));
4139566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Btype));
4149566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4159566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4169566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Btype));
417a6216909SToby Isaac     }
418a6216909SToby Isaac     if (D) {
419a6216909SToby Isaac       MPI_Datatype Dtype;
420a6216909SToby Isaac 
4219566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc * dim, origtype, &Dtype));
4229566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Dtype));
4239566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4249566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4259566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Dtype));
426a6216909SToby Isaac     }
427a6216909SToby Isaac     if (H) {
428a6216909SToby Isaac       MPI_Datatype Htype;
429a6216909SToby Isaac 
4309566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc * dim * dim, origtype, &Htype));
4319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Htype));
4329566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4339566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4349566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Htype));
435a6216909SToby Isaac     }
436a6216909SToby Isaac   }
4379566063dSJacob Faibussowitsch   PetscCall(PetscFree4(v, J, invJ, detJ));
4389566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBr, cellDr, cellHr));
4399566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBs, cellDs, cellHs));
4409566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellPoints, coordsReal, coordsRef));
4419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&pointType));
4429566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
443*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
444a6216909SToby Isaac }
445a6216909SToby Isaac 
446d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H)
447d71ae5a4SJacob Faibussowitsch {
448805e7170SToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
449805e7170SToby Isaac   PetscInt         h, imin;
450805e7170SToby Isaac   PetscInt         dim;
451805e7170SToby Isaac   PetscClassId     id;
452805e7170SToby Isaac   PetscQuadrature  quad = NULL;
453b7260050SToby Isaac   PetscInt         maxDegree;
454805e7170SToby Isaac   PetscFEGeom     *geom;
455805e7170SToby Isaac   PetscInt         Nq, Nc, dimC, qNc, N;
456805e7170SToby Isaac   PetscInt         numPoints;
457805e7170SToby Isaac   void            *qB = NULL, *qD = NULL, *qH = NULL;
458805e7170SToby Isaac   const PetscReal *weights;
459805e7170SToby Isaac   MPI_Datatype     mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL;
460805e7170SToby Isaac   PetscObject      disc;
461805e7170SToby Isaac   DMField          coordField;
462805e7170SToby Isaac 
463805e7170SToby Isaac   PetscFunctionBegin;
464805e7170SToby Isaac   Nc = field->numComponents;
4659566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dimC));
4669566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dim));
4679566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numPoints));
4689566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, NULL));
469805e7170SToby Isaac   for (h = 0; h < dsfield->height; h++) {
470805e7170SToby Isaac     PetscInt hEnd;
471805e7170SToby Isaac 
4729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
473805e7170SToby Isaac     if (imin < hEnd) break;
474805e7170SToby Isaac   }
475805e7170SToby Isaac   dim -= h;
4766858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
4779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
47808401ef6SPierre Jolivet   PetscCheck(id == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported");
4799566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(field->dm, &coordField));
4809566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
48148a46eb9SPierre Jolivet   if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
4829566063dSJacob Faibussowitsch   if (!quad) PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, &quad));
48328b400f6SJacob Faibussowitsch   PetscCheck(quad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages");
4849566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom));
4859566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights));
48608401ef6SPierre Jolivet   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components");
487805e7170SToby Isaac   N = numPoints * Nq * Nc;
4889566063dSJacob Faibussowitsch   if (B) PetscCall(DMGetWorkArray(field->dm, N, mpitype, &qB));
4899566063dSJacob Faibussowitsch   if (D) PetscCall(DMGetWorkArray(field->dm, N * dimC, mpitype, &qD));
4909566063dSJacob Faibussowitsch   if (H) PetscCall(DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
4919566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, pointIS, quad, type, qB, qD, qH));
492805e7170SToby Isaac   if (B) {
493805e7170SToby Isaac     PetscInt i, j, k;
494805e7170SToby Isaac 
495805e7170SToby Isaac     if (type == PETSC_SCALAR) {
496805e7170SToby Isaac       PetscScalar *sB  = (PetscScalar *)B;
497805e7170SToby Isaac       PetscScalar *sqB = (PetscScalar *)qB;
498805e7170SToby Isaac 
499805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
500805e7170SToby Isaac         PetscReal vol = 0.;
501805e7170SToby Isaac 
502ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) sB[i * Nc + j] = 0.;
503805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
504805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
505ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[(i * Nq + k) * Nc + j];
506805e7170SToby Isaac         }
50762bd480fSMatthew G. Knepley         for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol;
508805e7170SToby Isaac       }
509805e7170SToby Isaac     } else {
510805e7170SToby Isaac       PetscReal *rB  = (PetscReal *)B;
511805e7170SToby Isaac       PetscReal *rqB = (PetscReal *)qB;
512805e7170SToby Isaac 
513805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
514805e7170SToby Isaac         PetscReal vol = 0.;
515805e7170SToby Isaac 
516ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) rB[i * Nc + j] = 0.;
517805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
518805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
519ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) rB[i * Nc + j] += weights[k] * rqB[(i * Nq + k) * Nc + j];
520805e7170SToby Isaac         }
521805e7170SToby Isaac         for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol;
522805e7170SToby Isaac       }
523805e7170SToby Isaac     }
524805e7170SToby Isaac   }
525805e7170SToby Isaac   if (D) {
526805e7170SToby Isaac     PetscInt i, j, k, l, m;
527805e7170SToby Isaac 
528805e7170SToby Isaac     if (type == PETSC_SCALAR) {
529805e7170SToby Isaac       PetscScalar *sD  = (PetscScalar *)D;
530805e7170SToby Isaac       PetscScalar *sqD = (PetscScalar *)qD;
531805e7170SToby Isaac 
532805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
533805e7170SToby Isaac         PetscReal vol = 0.;
534805e7170SToby Isaac 
535ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) sD[i * Nc * dimC + j] = 0.;
536805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
537805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
538805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
539805e7170SToby Isaac             PetscScalar pD[3] = {0., 0., 0.};
540805e7170SToby Isaac 
541805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
542ad540459SPierre Jolivet               for (m = 0; m < dim; m++) pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * sqD[((i * Nq + k) * Nc + j) * dim + m];
543805e7170SToby Isaac             }
544ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
545805e7170SToby Isaac           }
546805e7170SToby Isaac         }
547805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol;
548805e7170SToby Isaac       }
549805e7170SToby Isaac     } else {
550805e7170SToby Isaac       PetscReal *rD  = (PetscReal *)D;
551805e7170SToby Isaac       PetscReal *rqD = (PetscReal *)qD;
552805e7170SToby Isaac 
553805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
554805e7170SToby Isaac         PetscReal vol = 0.;
555805e7170SToby Isaac 
556ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) rD[i * Nc * dimC + j] = 0.;
557805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
558805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
559805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
560805e7170SToby Isaac             PetscReal pD[3] = {0., 0., 0.};
561805e7170SToby Isaac 
562805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
563ad540459SPierre Jolivet               for (m = 0; m < dim; m++) pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * rqD[((i * Nq + k) * Nc + j) * dim + m];
564805e7170SToby Isaac             }
565ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
566805e7170SToby Isaac           }
567805e7170SToby Isaac         }
568805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol;
569805e7170SToby Isaac       }
570805e7170SToby Isaac     }
571805e7170SToby Isaac   }
572805e7170SToby Isaac   if (H) {
573805e7170SToby Isaac     PetscInt i, j, k, l, m, q, r;
574805e7170SToby Isaac 
575805e7170SToby Isaac     if (type == PETSC_SCALAR) {
576805e7170SToby Isaac       PetscScalar *sH  = (PetscScalar *)H;
577805e7170SToby Isaac       PetscScalar *sqH = (PetscScalar *)qH;
578805e7170SToby Isaac 
579805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
580805e7170SToby Isaac         PetscReal vol = 0.;
581805e7170SToby Isaac 
582ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) sH[i * Nc * dimC * dimC + j] = 0.;
583805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
5844d1a973fSToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
585805e7170SToby Isaac 
586805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
587805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
5889371c9d4SSatish Balay             PetscScalar pH[3][3] = {
5899371c9d4SSatish Balay               {0., 0., 0.},
5909371c9d4SSatish Balay               {0., 0., 0.},
5919371c9d4SSatish Balay               {0., 0., 0.}
5929371c9d4SSatish Balay             };
593805e7170SToby Isaac             const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC];
594805e7170SToby Isaac 
595805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
596805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
597805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
598ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r];
599805e7170SToby Isaac                 }
600805e7170SToby Isaac               }
601805e7170SToby Isaac             }
602805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
603ad540459SPierre Jolivet               for (m = 0; m < dimC; m++) sH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
604805e7170SToby Isaac             }
605805e7170SToby Isaac           }
606805e7170SToby Isaac         }
607805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol;
608805e7170SToby Isaac       }
609805e7170SToby Isaac     } else {
610805e7170SToby Isaac       PetscReal *rH  = (PetscReal *)H;
611805e7170SToby Isaac       PetscReal *rqH = (PetscReal *)qH;
612805e7170SToby Isaac 
613805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
614805e7170SToby Isaac         PetscReal vol = 0.;
615805e7170SToby Isaac 
616ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) rH[i * Nc * dimC * dimC + j] = 0.;
617805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
618805e7170SToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
619805e7170SToby Isaac 
620805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
621805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
6229371c9d4SSatish Balay             PetscReal pH[3][3] = {
6239371c9d4SSatish Balay               {0., 0., 0.},
6249371c9d4SSatish Balay               {0., 0., 0.},
6259371c9d4SSatish Balay               {0., 0., 0.}
6269371c9d4SSatish Balay             };
627805e7170SToby Isaac             const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC];
628805e7170SToby Isaac 
629805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
630805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
631805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
632ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r];
633805e7170SToby Isaac                 }
634805e7170SToby Isaac               }
635805e7170SToby Isaac             }
636805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
637ad540459SPierre Jolivet               for (m = 0; m < dimC; m++) rH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
638805e7170SToby Isaac             }
639805e7170SToby Isaac           }
640805e7170SToby Isaac         }
641805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol;
642805e7170SToby Isaac       }
643805e7170SToby Isaac     }
644805e7170SToby Isaac   }
6459566063dSJacob Faibussowitsch   if (B) PetscCall(DMRestoreWorkArray(field->dm, N, mpitype, &qB));
6469566063dSJacob Faibussowitsch   if (D) PetscCall(DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD));
6479566063dSJacob Faibussowitsch   if (H) PetscCall(DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
6489566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
6499566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
650*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
651805e7170SToby Isaac }
652805e7170SToby Isaac 
653d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
654d71ae5a4SJacob Faibussowitsch {
65551a454edSToby Isaac   DMField_DS  *dsfield;
65651a454edSToby Isaac   PetscObject  disc;
657741db25dSToby Isaac   PetscInt     h, imin, imax;
65851a454edSToby Isaac   PetscClassId id;
65951a454edSToby Isaac 
66051a454edSToby Isaac   PetscFunctionBegin;
66151a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
6629566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
663741db25dSToby Isaac   if (imin >= imax) {
66451aa0bf7SToby Isaac     h = 0;
66551aa0bf7SToby Isaac   } else {
666f99c8401SToby Isaac     for (h = 0; h < dsfield->height; h++) {
66751a454edSToby Isaac       PetscInt hEnd;
66851a454edSToby Isaac 
6699566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
67051a454edSToby Isaac       if (imin < hEnd) break;
67151a454edSToby Isaac     }
67251aa0bf7SToby Isaac   }
6736858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
6749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
67551a454edSToby Isaac   if (id == PETSCFE_CLASSID) {
67651a454edSToby Isaac     PetscFE    fe = (PetscFE)disc;
67751a454edSToby Isaac     PetscSpace sp;
67851a454edSToby Isaac 
6799566063dSJacob Faibussowitsch     PetscCall(PetscFEGetBasisSpace(fe, &sp));
6809566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetDegree(sp, minDegree, maxDegree));
68151a454edSToby Isaac   }
682*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68351a454edSToby Isaac }
68451a454edSToby Isaac 
685d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetFVQuadrature_Internal(DMField field, IS pointIS, PetscQuadrature *quad)
686d71ae5a4SJacob Faibussowitsch {
687a2aba86cSMatthew G. Knepley   DM              dm = field->dm;
688a2aba86cSMatthew G. Knepley   const PetscInt *points;
689a2aba86cSMatthew G. Knepley   DMPolytopeType  ct;
690a2aba86cSMatthew G. Knepley   PetscInt        dim, n;
691a2aba86cSMatthew G. Knepley   PetscBool       isplex;
692a2aba86cSMatthew G. Knepley 
693a2aba86cSMatthew G. Knepley   PetscFunctionBegin;
6949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex));
6959566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &n));
696a2aba86cSMatthew G. Knepley   if (isplex && n) {
6979566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
6989566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pointIS, &points));
6999566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, points[0], &ct));
700a2aba86cSMatthew G. Knepley     switch (ct) {
701a2aba86cSMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
702d71ae5a4SJacob Faibussowitsch     case DM_POLYTOPE_TETRAHEDRON:
703d71ae5a4SJacob Faibussowitsch       PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 1, -1.0, 1.0, quad));
704d71ae5a4SJacob Faibussowitsch       break;
705d71ae5a4SJacob Faibussowitsch     default:
706d71ae5a4SJacob Faibussowitsch       PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
707a2aba86cSMatthew G. Knepley     }
7089566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pointIS, &points));
7091baa6e33SBarry Smith   } else PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, quad));
710*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711a2aba86cSMatthew G. Knepley }
712a2aba86cSMatthew G. Knepley 
713d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad)
714d71ae5a4SJacob Faibussowitsch {
7151fdf1f07SMatthew G. Knepley   PetscInt     h, dim, imax, imin, cellHeight;
71651a454edSToby Isaac   DM           dm;
71751a454edSToby Isaac   DMField_DS  *dsfield;
71851a454edSToby Isaac   PetscObject  disc;
71951a454edSToby Isaac   PetscFE      fe;
72051a454edSToby Isaac   PetscClassId id;
72151a454edSToby Isaac 
72251a454edSToby Isaac   PetscFunctionBegin;
72351a454edSToby Isaac   dm      = field->dm;
72451a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
7259566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
7269566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
72751a454edSToby Isaac   for (h = 0; h <= dim; h++) {
72851a454edSToby Isaac     PetscInt hStart, hEnd;
72951a454edSToby Isaac 
7309566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
7316c041b70SSatish Balay     if (imax >= hStart && imin < hEnd) break;
73251a454edSToby Isaac   }
7339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
7341fdf1f07SMatthew G. Knepley   h -= cellHeight;
73551a454edSToby Isaac   *quad = NULL;
736f99c8401SToby Isaac   if (h < dsfield->height) {
7376858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
7389566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
739*3ba16761SJacob Faibussowitsch     if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
74051a454edSToby Isaac     fe = (PetscFE)disc;
7419566063dSJacob Faibussowitsch     PetscCall(PetscFEGetQuadrature(fe, quad));
7429566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)*quad));
74351a454edSToby Isaac   }
744*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74551a454edSToby Isaac }
74651a454edSToby Isaac 
747d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom)
748d71ae5a4SJacob Faibussowitsch {
7490145028aSToby Isaac   const PetscInt *points;
7500145028aSToby Isaac   PetscInt        p, dim, dE, numFaces, Nq;
751b7260050SToby Isaac   PetscInt        maxDegree;
7520145028aSToby Isaac   DMLabel         depthLabel;
7530145028aSToby Isaac   IS              cellIS;
7540145028aSToby Isaac   DM              dm = field->dm;
7550145028aSToby Isaac 
7560145028aSToby Isaac   PetscFunctionBegin;
7570145028aSToby Isaac   dim = geom->dim;
7580145028aSToby Isaac   dE  = geom->dimEmbed;
7599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
7609566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS));
7619566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(field, cellIS, NULL, &maxDegree));
7629566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointIS, &points));
7630145028aSToby Isaac   numFaces = geom->numCells;
7640145028aSToby Isaac   Nq       = geom->numPoints;
765f15274beSMatthew Knepley   /* First, set local faces and flip normals so that they are outward for the first supporting cell */
766f15274beSMatthew Knepley   for (p = 0; p < numFaces; p++) {
767f15274beSMatthew Knepley     PetscInt        point = points[p];
768f15274beSMatthew Knepley     PetscInt        suppSize, s, coneSize, c, numChildren;
769f15274beSMatthew Knepley     const PetscInt *supp, *cone, *ornt;
770f15274beSMatthew Knepley 
7719566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
77228b400f6SJacob Faibussowitsch     PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
7739566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, point, &suppSize));
77463a3b9bcSJacob Faibussowitsch     PetscCheck(suppSize <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, suppSize);
775f15274beSMatthew Knepley     if (!suppSize) continue;
7769566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, point, &supp));
777f15274beSMatthew Knepley     for (s = 0; s < suppSize; ++s) {
7789566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, supp[s], &coneSize));
7799566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, supp[s], &cone));
7809371c9d4SSatish Balay       for (c = 0; c < coneSize; ++c)
7819371c9d4SSatish Balay         if (cone[c] == point) break;
78263a3b9bcSJacob Faibussowitsch       PetscCheck(c != coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid connectivity: point %" PetscInt_FMT " not found in cone of support point %" PetscInt_FMT, point, supp[s]);
783f15274beSMatthew Knepley       geom->face[p][s] = c;
784f15274beSMatthew Knepley     }
7859566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeOrientation(dm, supp[0], &ornt));
786f15274beSMatthew Knepley     if (ornt[geom->face[p][0]] < 0) {
787f15274beSMatthew Knepley       PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d;
788f15274beSMatthew Knepley 
7899371c9d4SSatish Balay       for (q = 0; q < Np; ++q)
7909371c9d4SSatish Balay         for (d = 0; d < dE; ++d) geom->n[(p * Np + q) * dE + d] = -geom->n[(p * Np + q) * dE + d];
791f15274beSMatthew Knepley     }
792f15274beSMatthew Knepley   }
793b7260050SToby Isaac   if (maxDegree <= 1) {
7940145028aSToby Isaac     PetscInt     numCells, offset, *cells;
7950145028aSToby Isaac     PetscFEGeom *cellGeom;
7960145028aSToby Isaac     IS           suppIS;
7970145028aSToby Isaac 
7980145028aSToby Isaac     for (p = 0, numCells = 0; p < numFaces; p++) {
7990145028aSToby Isaac       PetscInt point = points[p];
8000145028aSToby Isaac       PetscInt numSupp, numChildren;
8010145028aSToby Isaac 
8029566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
80328b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
8049566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
80563a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
8060145028aSToby Isaac       numCells += numSupp;
8070145028aSToby Isaac     }
8089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numCells, &cells));
8090145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8100145028aSToby Isaac       PetscInt        point = points[p];
8110145028aSToby Isaac       PetscInt        numSupp, s;
8120145028aSToby Isaac       const PetscInt *supp;
8130145028aSToby Isaac 
8149566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
816ad540459SPierre Jolivet       for (s = 0; s < numSupp; s++, offset++) cells[offset] = supp[s];
8170145028aSToby Isaac     }
8189566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
8199566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateFEGeom(field, suppIS, quad, PETSC_FALSE, &cellGeom));
8200145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8210145028aSToby Isaac       PetscInt        point = points[p];
8220145028aSToby Isaac       PetscInt        numSupp, s, q;
8230145028aSToby Isaac       const PetscInt *supp;
8240145028aSToby Isaac 
8259566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8269566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
8270145028aSToby Isaac       for (s = 0; s < numSupp; s++, offset++) {
8280145028aSToby Isaac         for (q = 0; q < Nq * dE * dE; q++) {
82941418987SMatthew G. Knepley           geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
8300145028aSToby Isaac           geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
8310145028aSToby Isaac         }
83241418987SMatthew G. Knepley         for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
8330145028aSToby Isaac       }
8340145028aSToby Isaac     }
8359566063dSJacob Faibussowitsch     PetscCall(PetscFEGeomDestroy(&cellGeom));
8369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&suppIS));
8379566063dSJacob Faibussowitsch     PetscCall(PetscFree(cells));
8380145028aSToby Isaac   } else {
8396858538eSMatthew G. Knepley     DMField_DS    *dsfield = (DMField_DS *)field->data;
8400145028aSToby Isaac     PetscObject    faceDisc, cellDisc;
8410145028aSToby Isaac     PetscClassId   faceId, cellId;
8420145028aSToby Isaac     PetscDualSpace dsp;
8430145028aSToby Isaac     DM             K;
844ba2698f1SMatthew G. Knepley     DMPolytopeType ct;
8450145028aSToby Isaac     PetscInt(*co)[2][3];
8460145028aSToby Isaac     PetscInt        coneSize;
8470145028aSToby Isaac     PetscInt      **counts;
8480145028aSToby Isaac     PetscInt        f, i, o, q, s;
8490145028aSToby Isaac     const PetscInt *coneK;
850ba2698f1SMatthew G. Knepley     PetscInt        eStart, minOrient, maxOrient, numOrient;
8510145028aSToby Isaac     PetscInt       *orients;
8520145028aSToby Isaac     PetscReal     **orientPoints;
8530145028aSToby Isaac     PetscReal      *cellPoints;
8540145028aSToby Isaac     PetscReal      *dummyWeights;
8550145028aSToby Isaac     PetscQuadrature cellQuad = NULL;
8560145028aSToby Isaac 
8576858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 1, dsfield->disc, &faceDisc));
8586858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
8599566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(faceDisc, &faceId));
8609566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(cellDisc, &cellId));
8611dca8a05SBarry Smith     PetscCheck(faceId == PETSCFE_CLASSID && cellId == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported");
8629566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace((PetscFE)cellDisc, &dsp));
8639566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(dsp, &K));
8649566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(K, 1, &eStart, NULL));
8659566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(K, eStart, &ct));
8669566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(K, 0, &coneSize));
8679566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
8689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numFaces, &co, coneSize, &counts));
8699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dE * Nq, &cellPoints));
8709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nq, &dummyWeights));
8719566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad));
8729566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights));
8730145028aSToby Isaac     minOrient = PETSC_MAX_INT;
8740145028aSToby Isaac     maxOrient = PETSC_MIN_INT;
8750145028aSToby Isaac     for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */
8760145028aSToby Isaac       PetscInt        point = points[p];
8770145028aSToby Isaac       PetscInt        numSupp, numChildren;
8780145028aSToby Isaac       const PetscInt *supp;
8790145028aSToby Isaac 
8809566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
88128b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
8829566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
88363a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
8849566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
8850145028aSToby Isaac       for (s = 0; s < numSupp; s++) {
8860145028aSToby Isaac         PetscInt        cell = supp[s];
8870145028aSToby Isaac         PetscInt        numCone;
8880145028aSToby Isaac         const PetscInt *cone, *orient;
8890145028aSToby Isaac 
8909566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(dm, cell, &numCone));
89108401ef6SPierre Jolivet         PetscCheck(numCone == coneSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element");
8929566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(dm, cell, &cone));
8939566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeOrientation(dm, cell, &orient));
8940145028aSToby Isaac         for (f = 0; f < coneSize; f++) {
8950145028aSToby Isaac           if (cone[f] == point) break;
8960145028aSToby Isaac         }
8970145028aSToby Isaac         co[p][s][0] = f;
8980145028aSToby Isaac         co[p][s][1] = orient[f];
8990145028aSToby Isaac         co[p][s][2] = cell;
9000145028aSToby Isaac         minOrient   = PetscMin(minOrient, orient[f]);
901cd93b0e1SMatthew G. Knepley         maxOrient   = PetscMax(maxOrient, orient[f]);
9020145028aSToby Isaac       }
9030145028aSToby Isaac       for (; s < 2; s++) {
9040145028aSToby Isaac         co[p][s][0] = -1;
9050145028aSToby Isaac         co[p][s][1] = -1;
9060145028aSToby Isaac         co[p][s][2] = -1;
9070145028aSToby Isaac       }
9080145028aSToby Isaac     }
9090145028aSToby Isaac     numOrient = maxOrient + 1 - minOrient;
9109566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
9110145028aSToby Isaac     /* count all (face,orientation) doubles that appear */
9129566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(numOrient, &orients, numOrient, &orientPoints));
9139566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscCalloc1(numOrient + 1, &counts[f]));
9140145028aSToby Isaac     for (p = 0; p < numFaces; p++) {
9150145028aSToby Isaac       for (s = 0; s < 2; s++) {
9160145028aSToby Isaac         if (co[p][s][0] >= 0) {
9170145028aSToby Isaac           counts[co[p][s][0]][co[p][s][1] - minOrient]++;
9180145028aSToby Isaac           orients[co[p][s][1] - minOrient]++;
9190145028aSToby Isaac         }
9200145028aSToby Isaac       }
9210145028aSToby Isaac     }
9220145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
9230145028aSToby Isaac       if (orients[o]) {
9240145028aSToby Isaac         PetscInt orient = o + minOrient;
9250145028aSToby Isaac         PetscInt q;
9260145028aSToby Isaac 
9279566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Nq * dim, &orientPoints[o]));
9280145028aSToby Isaac         /* rotate the quadrature points appropriately */
929ba2698f1SMatthew G. Knepley         switch (ct) {
930d71ae5a4SJacob Faibussowitsch         case DM_POLYTOPE_POINT:
931d71ae5a4SJacob Faibussowitsch           break;
932ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_SEGMENT:
9330145028aSToby Isaac           if (orient == -2 || orient == 1) {
934ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = -geom->xi[q];
9350145028aSToby Isaac           } else {
936ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = geom->xi[q];
9370145028aSToby Isaac           }
9380145028aSToby Isaac           break;
939ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:
9400145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9410145028aSToby Isaac             PetscReal lambda[3];
9420145028aSToby Isaac             PetscReal lambdao[3];
9430145028aSToby Isaac 
9440145028aSToby Isaac             /* convert to barycentric */
9450145028aSToby Isaac             lambda[0] = -(geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.;
9460145028aSToby Isaac             lambda[1] = (geom->xi[2 * q] + 1.) / 2.;
9470145028aSToby Isaac             lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.;
9480145028aSToby Isaac             if (orient >= 0) {
949ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(orient + i) % 3];
9500145028aSToby Isaac             } else {
951ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(-(orient + i) + 3) % 3];
9520145028aSToby Isaac             }
9530145028aSToby Isaac             /* convert to coordinates */
9540145028aSToby Isaac             orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1];
9550145028aSToby Isaac             orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2];
9560145028aSToby Isaac           }
9570145028aSToby Isaac           break;
958ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL:
9590145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9600145028aSToby Isaac             PetscReal xi[2], xio[2];
9610145028aSToby Isaac             PetscInt  oabs = (orient >= 0) ? orient : -(orient + 1);
9620145028aSToby Isaac 
9630145028aSToby Isaac             xi[0] = geom->xi[2 * q];
9640145028aSToby Isaac             xi[1] = geom->xi[2 * q + 1];
9650145028aSToby Isaac             switch (oabs) {
9660145028aSToby Isaac             case 1:
9670145028aSToby Isaac               xio[0] = xi[1];
9680145028aSToby Isaac               xio[1] = -xi[0];
9690145028aSToby Isaac               break;
970f4d061e9SPierre Jolivet             case 2:
971f4d061e9SPierre Jolivet               xio[0] = -xi[0];
972f4d061e9SPierre Jolivet               xio[1] = -xi[1];
973f4d061e9SPierre Jolivet               break;
974f4d061e9SPierre Jolivet             case 3:
975f4d061e9SPierre Jolivet               xio[0] = -xi[1];
976f4d061e9SPierre Jolivet               xio[1] = xi[0];
977f4d061e9SPierre Jolivet               break;
978aa1371cbSToby Isaac             case 0:
979aa1371cbSToby Isaac             default:
980aa1371cbSToby Isaac               xio[0] = xi[0];
981aa1371cbSToby Isaac               xio[1] = xi[1];
982aa1371cbSToby Isaac               break;
9830145028aSToby Isaac             }
984ad540459SPierre Jolivet             if (orient < 0) xio[0] = -xio[0];
9850145028aSToby Isaac             orientPoints[o][2 * q + 0] = xio[0];
9860145028aSToby Isaac             orientPoints[o][2 * q + 1] = xio[1];
9870145028aSToby Isaac           }
9880145028aSToby Isaac           break;
989d71ae5a4SJacob Faibussowitsch         default:
990d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s not yet supported", DMPolytopeTypes[ct]);
9910145028aSToby Isaac         }
9920145028aSToby Isaac       }
9930145028aSToby Isaac     }
9940145028aSToby Isaac     for (f = 0; f < coneSize; f++) {
9950145028aSToby Isaac       PetscInt  face = coneK[f];
9964d1a973fSToby Isaac       PetscReal v0[3];
9972d89661fSMatthew G. Knepley       PetscReal J[9], detJ;
9980145028aSToby Isaac       PetscInt  numCells, offset;
9990145028aSToby Isaac       PetscInt *cells;
10000145028aSToby Isaac       IS        suppIS;
10010145028aSToby Isaac 
10029566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ));
10030145028aSToby Isaac       for (o = 0; o <= numOrient; o++) {
10040145028aSToby Isaac         PetscFEGeom *cellGeom;
10050145028aSToby Isaac 
10060145028aSToby Isaac         if (!counts[f][o]) continue;
10070145028aSToby Isaac         /* If this (face,orientation) double appears,
10080145028aSToby Isaac          * convert the face quadrature points into volume quadrature points */
10090145028aSToby Isaac         for (q = 0; q < Nq; q++) {
10100145028aSToby Isaac           PetscReal xi0[3] = {-1., -1., -1.};
10110145028aSToby Isaac 
10120145028aSToby Isaac           CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]);
10130145028aSToby Isaac         }
10140145028aSToby Isaac         for (p = 0, numCells = 0; p < numFaces; p++) {
10150145028aSToby Isaac           for (s = 0; s < 2; s++) {
10160145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++;
10170145028aSToby Isaac           }
10180145028aSToby Isaac         }
10199566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(numCells, &cells));
10200145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10210145028aSToby Isaac           for (s = 0; s < 2; s++) {
1022ad540459SPierre Jolivet             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) cells[offset++] = co[p][s][2];
10230145028aSToby Isaac           }
10240145028aSToby Isaac         }
10259566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
10269566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FALSE, &cellGeom));
10270145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10280145028aSToby Isaac           for (s = 0; s < 2; s++) {
10290145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
10300145028aSToby Isaac               for (q = 0; q < Nq * dE * dE; q++) {
103141418987SMatthew G. Knepley                 geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
10320145028aSToby Isaac                 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
10330145028aSToby Isaac               }
103441418987SMatthew G. Knepley               for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
10350145028aSToby Isaac               offset++;
10360145028aSToby Isaac             }
10370145028aSToby Isaac           }
10380145028aSToby Isaac         }
10399566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&cellGeom));
10409566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&suppIS));
10419566063dSJacob Faibussowitsch         PetscCall(PetscFree(cells));
10420145028aSToby Isaac       }
10430145028aSToby Isaac     }
10440145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
104548a46eb9SPierre Jolivet       if (orients[o]) PetscCall(PetscFree(orientPoints[o]));
10460145028aSToby Isaac     }
10479566063dSJacob Faibussowitsch     PetscCall(PetscFree2(orients, orientPoints));
10489566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&cellQuad));
10499566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscFree(counts[f]));
10509566063dSJacob Faibussowitsch     PetscCall(PetscFree2(co, counts));
10510145028aSToby Isaac   }
10529566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointIS, &points));
10539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
1054*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10550145028aSToby Isaac }
10560145028aSToby Isaac 
1057d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldInitialize_DS(DMField field)
1058d71ae5a4SJacob Faibussowitsch {
105951a454edSToby Isaac   PetscFunctionBegin;
106051a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DS;
106151a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DS;
106251a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DS;
1063805e7170SToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DS;
1064b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DS;
106551a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS;
106651a454edSToby Isaac   field->ops->view                    = DMFieldView_DS;
10670145028aSToby Isaac   field->ops->computeFaceData         = DMFieldComputeFaceData_DS;
1068*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106951a454edSToby Isaac }
107051a454edSToby Isaac 
1071d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field)
1072d71ae5a4SJacob Faibussowitsch {
107351a454edSToby Isaac   DMField_DS *dsfield;
107451a454edSToby Isaac 
107551a454edSToby Isaac   PetscFunctionBegin;
10764dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dsfield));
107751a454edSToby Isaac   field->data = dsfield;
10789566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_DS(field));
1079*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108051a454edSToby Isaac }
108151a454edSToby Isaac 
1082d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDSWithDG(DM dm, DM dmDG, PetscInt fieldNum, Vec vec, Vec vecDG, DMField *field)
1083d71ae5a4SJacob Faibussowitsch {
108451a454edSToby Isaac   DMField      b;
108551a454edSToby Isaac   DMField_DS  *dsfield;
10866858538eSMatthew G. Knepley   PetscObject  disc = NULL, discDG = NULL;
10876858538eSMatthew G. Knepley   PetscSection section;
108851a454edSToby Isaac   PetscBool    isContainer   = PETSC_FALSE;
108951a454edSToby Isaac   PetscClassId id            = -1;
109051a454edSToby Isaac   PetscInt     numComponents = -1, dsNumFields;
109151a454edSToby Isaac 
109251a454edSToby Isaac   PetscFunctionBegin;
10939566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
10949566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetFieldComponents(section, fieldNum, &numComponents));
10959566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &dsNumFields));
10969566063dSJacob Faibussowitsch   if (dsNumFields) PetscCall(DMGetField(dm, fieldNum, NULL, &disc));
10976858538eSMatthew G. Knepley   if (dsNumFields && dmDG) {
10986858538eSMatthew G. Knepley     PetscCall(DMGetField(dmDG, fieldNum, NULL, &discDG));
10996858538eSMatthew G. Knepley     PetscCall(PetscObjectReference(discDG));
11006858538eSMatthew G. Knepley   }
110151a454edSToby Isaac   if (disc) {
11029566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
110351a454edSToby Isaac     isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE;
110451a454edSToby Isaac   }
110551a454edSToby Isaac   if (!disc || isContainer) {
110651a454edSToby Isaac     MPI_Comm       comm = PetscObjectComm((PetscObject)dm);
110751a454edSToby Isaac     PetscFE        fe;
11082df84da0SMatthew G. Knepley     DMPolytopeType ct, locct = DM_POLYTOPE_UNKNOWN;
11092df84da0SMatthew G. Knepley     PetscInt       dim, cStart, cEnd, cellHeight;
111051a454edSToby Isaac 
11119566063dSJacob Faibussowitsch     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
11129566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
11139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
11149566063dSJacob Faibussowitsch     if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &locct));
11159566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&locct, &ct, 1, MPI_INT, MPI_MIN, comm));
11169566063dSJacob Faibussowitsch     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, numComponents, ct, 1, PETSC_DETERMINE, &fe));
11179566063dSJacob Faibussowitsch     PetscCall(PetscFEViewFromOptions(fe, NULL, "-field_fe_view"));
111851a454edSToby Isaac     disc = (PetscObject)fe;
11191baa6e33SBarry Smith   } else PetscCall(PetscObjectReference(disc));
11209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
11216858538eSMatthew G. Knepley   if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)disc, &numComponents));
11226858538eSMatthew G. Knepley   else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine number of discretization components");
11239566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm, numComponents, DMFIELD_VERTEX, &b));
11249566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b, DMFIELDDS));
112551a454edSToby Isaac   dsfield           = (DMField_DS *)b->data;
112651a454edSToby Isaac   dsfield->fieldNum = fieldNum;
11279566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dsfield->height));
1128f99c8401SToby Isaac   dsfield->height++;
11299566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(dsfield->height, &dsfield->disc));
113051a454edSToby Isaac   dsfield->disc[0] = disc;
11319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)vec));
113251a454edSToby Isaac   dsfield->vec = vec;
11336858538eSMatthew G. Knepley   if (dmDG) {
11346858538eSMatthew G. Knepley     dsfield->dmDG = dmDG;
11356858538eSMatthew G. Knepley     PetscCall(PetscCalloc1(dsfield->height, &dsfield->discDG));
11366858538eSMatthew G. Knepley     dsfield->discDG[0] = discDG;
11376858538eSMatthew G. Knepley     PetscCall(PetscObjectReference((PetscObject)vecDG));
11386858538eSMatthew G. Knepley     dsfield->vecDG = vecDG;
11396858538eSMatthew G. Knepley   }
114051a454edSToby Isaac   *field = b;
1141*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11426858538eSMatthew G. Knepley }
114351a454edSToby Isaac 
1144d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec, DMField *field)
1145d71ae5a4SJacob Faibussowitsch {
11466858538eSMatthew G. Knepley   PetscFunctionBegin;
11476858538eSMatthew G. Knepley   PetscCall(DMFieldCreateDSWithDG(dm, NULL, fieldNum, vec, NULL, field));
1148*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114951a454edSToby Isaac }
1150