xref: /petsc/src/dm/field/impls/ds/dmfieldds.c (revision 6a210b70a61472cf8733db59e8a3fa68f6578c01)
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));
323ba16761SJacob 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));
543ba16761SJacob 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;
60313f16f2SMatthew G. Knepley   PetscCheck(height >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Height %" PetscInt_FMT " must be non-negative", height);
616858538eSMatthew G. Knepley   if (!discList[height]) {
6251a454edSToby Isaac     PetscClassId id;
6351a454edSToby Isaac 
646858538eSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(discList[0], &id));
656858538eSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateHeightTrace((PetscFE)discList[0], height, (PetscFE *)&discList[height]));
6651a454edSToby Isaac   }
676858538eSMatthew G. Knepley   *disc = discList[height];
683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6951a454edSToby Isaac }
7051a454edSToby Isaac 
7162bd480fSMatthew G. Knepley /* y[m,c] = A[m,n,c] . b[n] */
7251a454edSToby Isaac #define DMFieldDSdot(y, A, b, m, n, c, cast) \
7351a454edSToby Isaac   do { \
7451a454edSToby Isaac     PetscInt _i, _j, _k; \
7551a454edSToby Isaac     for (_i = 0; _i < (m); _i++) { \
76ad540459SPierre Jolivet       for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] = 0.; \
7751a454edSToby Isaac       for (_j = 0; _j < (n); _j++) { \
78ad540459SPierre Jolivet         for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \
7951a454edSToby Isaac       } \
8051a454edSToby Isaac     } \
8151a454edSToby Isaac   } while (0)
8251a454edSToby Isaac 
836858538eSMatthew G. Knepley /*
846858538eSMatthew 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().
856858538eSMatthew G. Knepley */
8666976f2fSJacob Faibussowitsch static PetscErrorCode DMFieldGetClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[])
87d71ae5a4SJacob Faibussowitsch {
886858538eSMatthew G. Knepley   DMField_DS        *dsfield = (DMField_DS *)field->data;
896858538eSMatthew G. Knepley   DM                 fdm     = dsfield->dmDG;
906858538eSMatthew G. Knepley   PetscSection       s       = NULL;
916858538eSMatthew G. Knepley   const PetscScalar *cvalues;
926858538eSMatthew G. Knepley   PetscInt           pStart, pEnd;
936858538eSMatthew G. Knepley 
946858538eSMatthew G. Knepley   PetscFunctionBeginHot;
956858538eSMatthew G. Knepley   *isDG   = PETSC_FALSE;
966858538eSMatthew G. Knepley   *Nc     = 0;
976858538eSMatthew G. Knepley   *array  = NULL;
986858538eSMatthew G. Knepley   *values = NULL;
996858538eSMatthew G. Knepley   /* Check for cellwise section */
1006858538eSMatthew G. Knepley   if (fdm) PetscCall(DMGetLocalSection(fdm, &s));
1016858538eSMatthew G. Knepley   if (!s) goto cg;
1026858538eSMatthew G. Knepley   /* Check that the cell exists in the cellwise section */
1036858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
1046858538eSMatthew G. Knepley   if (cell < pStart || cell >= pEnd) goto cg;
1056858538eSMatthew G. Knepley   /* Check for cellwise coordinates for this cell */
1066858538eSMatthew G. Knepley   PetscCall(PetscSectionGetDof(s, cell, Nc));
1076858538eSMatthew G. Knepley   if (!*Nc) goto cg;
1086858538eSMatthew G. Knepley   /* Check for cellwise coordinates */
1096858538eSMatthew G. Knepley   if (!dsfield->vecDG) goto cg;
1106858538eSMatthew G. Knepley   /* Get cellwise coordinates */
1116858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(dsfield->vecDG, array));
1126858538eSMatthew G. Knepley   PetscCall(DMPlexPointLocalRead(fdm, cell, *array, &cvalues));
1136858538eSMatthew G. Knepley   PetscCall(DMGetWorkArray(fdm, *Nc, MPIU_SCALAR, values));
1146858538eSMatthew G. Knepley   PetscCall(PetscArraycpy(*values, cvalues, *Nc));
1156858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(dsfield->vecDG, array));
1166858538eSMatthew G. Knepley   *isDG = PETSC_TRUE;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1186858538eSMatthew G. Knepley cg:
1196858538eSMatthew G. Knepley   /* Use continuous values */
1206858538eSMatthew G. Knepley   PetscCall(DMFieldGetDM(field, &fdm));
1216858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(fdm, &s));
1226267c18cSMatthew G. Knepley   PetscCall(PetscSectionGetField(s, dsfield->fieldNum, &s));
1236858538eSMatthew G. Knepley   PetscCall(DMPlexVecGetClosure(fdm, s, dsfield->vec, cell, Nc, values));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1256858538eSMatthew G. Knepley }
1266858538eSMatthew G. Knepley 
12766976f2fSJacob Faibussowitsch static PetscErrorCode DMFieldRestoreClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[])
128d71ae5a4SJacob Faibussowitsch {
1296858538eSMatthew G. Knepley   DMField_DS  *dsfield = (DMField_DS *)field->data;
1306858538eSMatthew G. Knepley   DM           fdm;
1316858538eSMatthew G. Knepley   PetscSection s;
1326858538eSMatthew G. Knepley 
1336858538eSMatthew G. Knepley   PetscFunctionBeginHot;
1346858538eSMatthew G. Knepley   if (*isDG) {
1356858538eSMatthew G. Knepley     PetscCall(DMRestoreWorkArray(dsfield->dmDG, *Nc, MPIU_SCALAR, values));
1366858538eSMatthew G. Knepley   } else {
1376858538eSMatthew G. Knepley     PetscCall(DMFieldGetDM(field, &fdm));
1386858538eSMatthew G. Knepley     PetscCall(DMGetLocalSection(fdm, &s));
1396858538eSMatthew G. Knepley     PetscCall(DMPlexVecRestoreClosure(fdm, s, dsfield->vec, cell, Nc, (PetscScalar **)values));
1406858538eSMatthew G. Knepley   }
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1426858538eSMatthew G. Knepley }
1436858538eSMatthew G. Knepley 
144ef0bb6c7SMatthew G. Knepley /* TODO: Reorganize interface so that I can reuse a tabulation rather than mallocing each time */
145d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H)
146d71ae5a4SJacob Faibussowitsch {
14751a454edSToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
14851a454edSToby Isaac   DM               dm;
14951a454edSToby Isaac   PetscObject      disc;
15051a454edSToby Isaac   PetscClassId     classid;
15151a454edSToby Isaac   PetscInt         nq, nc, dim, meshDim, numCells;
15251a454edSToby Isaac   PetscSection     section;
15351a454edSToby Isaac   const PetscReal *qpoints;
15451a454edSToby Isaac   PetscBool        isStride;
15551a454edSToby Isaac   const PetscInt  *points = NULL;
15651a454edSToby Isaac   PetscInt         sfirst = -1, stride = -1;
15751a454edSToby Isaac 
15851a454edSToby Isaac   PetscFunctionBeginHot;
15951a454edSToby Isaac   dm = field->dm;
16051a454edSToby Isaac   nc = field->numComponents;
1619566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &nq, &qpoints, NULL));
1626858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, dsfield->height - 1 - dim, dsfield->disc, &disc));
1639566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &meshDim));
1649566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1659566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetField(section, dsfield->fieldNum, &section));
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &classid));
16751a454edSToby Isaac   /* TODO: batch */
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
1699566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
1701baa6e33SBarry Smith   if (isStride) PetscCall(ISStrideGetInfo(pointIS, &sfirst, &stride));
1711baa6e33SBarry Smith   else PetscCall(ISGetIndices(pointIS, &points));
17251a454edSToby Isaac   if (classid == PETSCFE_CLASSID) {
17351a454edSToby Isaac     PetscFE         fe = (PetscFE)disc;
17451a454edSToby Isaac     PetscInt        feDim, i;
175ef0bb6c7SMatthew G. Knepley     PetscInt        K = H ? 2 : (D ? 1 : (B ? 0 : -1));
176ef0bb6c7SMatthew G. Knepley     PetscTabulation T;
17751a454edSToby Isaac 
1789566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDimension(fe, &feDim));
1799566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fe, 1, nq, qpoints, K, &T));
18051a454edSToby Isaac     for (i = 0; i < numCells; i++) {
18151a454edSToby Isaac       PetscInt           c = isStride ? (sfirst + i * stride) : points[i];
1822710589cSToby Isaac       PetscInt           closureSize;
1836858538eSMatthew G. Knepley       const PetscScalar *array;
1842710589cSToby Isaac       PetscScalar       *elem = NULL;
1856858538eSMatthew G. Knepley       PetscBool          isDG;
18651a454edSToby Isaac 
1876858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
18851a454edSToby Isaac       if (B) {
18962bd480fSMatthew G. Knepley         /* field[c] = T[q,b,c] . coef[b], so v[c] = T[q,b,c] . coords[b] */
19051a454edSToby Isaac         if (type == PETSC_SCALAR) {
19151a454edSToby Isaac           PetscScalar *cB = &((PetscScalar *)B)[nc * nq * i];
19251a454edSToby Isaac 
193ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
19451a454edSToby Isaac         } else {
19551a454edSToby Isaac           PetscReal *cB = &((PetscReal *)B)[nc * nq * i];
19651a454edSToby Isaac 
197ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
19851a454edSToby Isaac         }
19951a454edSToby Isaac       }
20051a454edSToby Isaac       if (D) {
20151a454edSToby Isaac         if (type == PETSC_SCALAR) {
20251a454edSToby Isaac           PetscScalar *cD = &((PetscScalar *)D)[nc * nq * dim * i];
20351a454edSToby Isaac 
204ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
20551a454edSToby Isaac         } else {
20651a454edSToby Isaac           PetscReal *cD = &((PetscReal *)D)[nc * nq * dim * i];
20751a454edSToby Isaac 
208ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
20951a454edSToby Isaac         }
21051a454edSToby Isaac       }
21151a454edSToby Isaac       if (H) {
21251a454edSToby Isaac         if (type == PETSC_SCALAR) {
21351a454edSToby Isaac           PetscScalar *cH = &((PetscScalar *)H)[nc * nq * dim * dim * i];
21451a454edSToby Isaac 
215ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
21651a454edSToby Isaac         } else {
21751a454edSToby Isaac           PetscReal *cH = &((PetscReal *)H)[nc * nq * dim * dim * i];
21851a454edSToby Isaac 
219ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
22051a454edSToby Isaac         }
22151a454edSToby Isaac       }
2226858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
22351a454edSToby Isaac     }
2249566063dSJacob Faibussowitsch     PetscCall(PetscTabulationDestroy(&T));
2252da392ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented");
22648a46eb9SPierre Jolivet   if (!isStride) PetscCall(ISRestoreIndices(pointIS, &points));
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22851a454edSToby Isaac }
22951a454edSToby Isaac 
230d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
231d71ae5a4SJacob Faibussowitsch {
232a6216909SToby Isaac   DMField_DS        *dsfield = (DMField_DS *)field->data;
233a6216909SToby Isaac   PetscSF            cellSF  = NULL;
234a6216909SToby Isaac   const PetscSFNode *cells;
235a6216909SToby Isaac   PetscInt           c, nFound, numCells, feDim, nc;
236a6216909SToby Isaac   const PetscInt    *cellDegrees;
237a6216909SToby Isaac   const PetscScalar *pointsArray;
238a6216909SToby Isaac   PetscScalar       *cellPoints;
239a6216909SToby Isaac   PetscInt           gatherSize, gatherMax;
240a6216909SToby Isaac   PetscInt           dim, dimR, offset;
241a6216909SToby Isaac   MPI_Datatype       pointType;
242a6216909SToby Isaac   PetscObject        cellDisc;
243a6216909SToby Isaac   PetscFE            cellFE;
244a6216909SToby Isaac   PetscClassId       discID;
245a6216909SToby Isaac   PetscReal         *coordsReal, *coordsRef;
246a6216909SToby Isaac   PetscSection       section;
247a6216909SToby Isaac   PetscScalar       *cellBs = NULL, *cellDs = NULL, *cellHs = NULL;
248a6216909SToby Isaac   PetscReal         *cellBr = NULL, *cellDr = NULL, *cellHr = NULL;
2494cbe5319SToby Isaac   PetscReal         *v, *J, *invJ, *detJ;
250a6216909SToby Isaac 
251a6216909SToby Isaac   PetscFunctionBegin;
252a6216909SToby Isaac   nc = field->numComponents;
2539566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(field->dm, &section));
2546858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
2559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(cellDisc, &discID));
25608401ef6SPierre Jolivet   PetscCheck(discID == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization type not supported");
257a6216909SToby Isaac   cellFE = (PetscFE)cellDisc;
2589566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDimension(cellFE, &feDim));
2599566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dim));
2609566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dimR));
2619566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF));
2629566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells));
263ad540459SPierre 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);
2649566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(cellSF, &cellDegrees));
2659566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(cellSF, &cellDegrees));
266a6216909SToby Isaac   for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) {
267a6216909SToby Isaac     gatherMax = PetscMax(gatherMax, cellDegrees[c]);
268a6216909SToby Isaac     gatherSize += cellDegrees[c];
269a6216909SToby Isaac   }
2709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(gatherSize * dim, &cellPoints, gatherMax * dim, &coordsReal, gatherMax * dimR, &coordsRef));
2719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(gatherMax * dimR, &v, gatherMax * dimR * dimR, &J, gatherMax * dimR * dimR, &invJ, gatherMax, &detJ));
2721baa6e33SBarry 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));
2731baa6e33SBarry Smith   else PetscCall(PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr));
274a6216909SToby Isaac 
2756497c311SBarry Smith   PetscCallMPI(MPI_Type_contiguous((PetscMPIInt)dim, MPIU_SCALAR, &pointType));
2769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&pointType));
2779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points, &pointsArray));
2789566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints));
2799566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints));
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points, &pointsArray));
281a6216909SToby Isaac   for (c = 0, offset = 0; c < numCells; c++) {
282a6216909SToby Isaac     PetscInt nq = cellDegrees[c], p;
283a6216909SToby Isaac 
284a6216909SToby Isaac     if (nq) {
285ef0bb6c7SMatthew G. Knepley       PetscInt           K = H ? 2 : (D ? 1 : (B ? 0 : -1));
286ef0bb6c7SMatthew G. Knepley       PetscTabulation    T;
2876858538eSMatthew G. Knepley       PetscQuadrature    quad;
2886858538eSMatthew G. Knepley       const PetscScalar *array;
289a6216909SToby Isaac       PetscScalar       *elem = NULL;
2904cbe5319SToby Isaac       PetscReal         *quadPoints;
2916858538eSMatthew G. Knepley       PetscBool          isDG;
2926858538eSMatthew G. Knepley       PetscInt           closureSize, d, e, f, g;
293a6216909SToby Isaac 
2944d1a973fSToby Isaac       for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]);
2959566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef));
2969566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(cellFE, 1, nq, coordsRef, K, &T));
2979566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad));
2989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dimR * nq, &quadPoints));
2994cbe5319SToby Isaac       for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p];
3009566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL));
3019566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ));
3029566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
3036858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
304a6216909SToby Isaac       if (B) {
305a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
306a6216909SToby Isaac           PetscScalar *cB = &cellBs[nc * offset];
307a6216909SToby Isaac 
308ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
309a6216909SToby Isaac         } else {
310a6216909SToby Isaac           PetscReal *cB = &cellBr[nc * offset];
311a6216909SToby Isaac 
312ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
313a6216909SToby Isaac         }
314a6216909SToby Isaac       }
315a6216909SToby Isaac       if (D) {
316a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
317a6216909SToby Isaac           PetscScalar *cD = &cellDs[nc * dim * offset];
318a6216909SToby Isaac 
319ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
3204cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3214cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3224d1a973fSToby Isaac               PetscScalar vs[3];
3234d1a973fSToby Isaac 
3244cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3254d1a973fSToby Isaac                 vs[d] = 0.;
326ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3274cbe5319SToby Isaac               }
328ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = vs[d];
3294cbe5319SToby Isaac             }
3304cbe5319SToby Isaac           }
331a6216909SToby Isaac         } else {
332a6216909SToby Isaac           PetscReal *cD = &cellDr[nc * dim * offset];
333a6216909SToby Isaac 
334ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
3354cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3364cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3374cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3384cbe5319SToby Isaac                 v[d] = 0.;
339ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3404cbe5319SToby Isaac               }
341ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = v[d];
3424cbe5319SToby Isaac             }
3434cbe5319SToby Isaac           }
344a6216909SToby Isaac         }
345a6216909SToby Isaac       }
346a6216909SToby Isaac       if (H) {
347a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
348a6216909SToby Isaac           PetscScalar *cH = &cellHs[nc * dim * dim * offset];
349a6216909SToby Isaac 
350ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
3514cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3524cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3534d1a973fSToby Isaac               PetscScalar vs[3];
3544d1a973fSToby Isaac 
3554cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3564d1a973fSToby Isaac                 vs[d] = 0.;
357ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3584cbe5319SToby Isaac               }
359ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = vs[d];
3604cbe5319SToby Isaac             }
3614cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3624cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3634d1a973fSToby Isaac                 PetscScalar vs[3];
3644d1a973fSToby Isaac 
3654cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3664d1a973fSToby Isaac                   vs[d] = 0.;
367ad540459SPierre Jolivet                   for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3684cbe5319SToby Isaac                 }
369ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d];
3704cbe5319SToby Isaac               }
3714cbe5319SToby Isaac             }
3724cbe5319SToby Isaac           }
373a6216909SToby Isaac         } else {
374a6216909SToby Isaac           PetscReal *cH = &cellHr[nc * dim * dim * offset];
375a6216909SToby Isaac 
376ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
3774cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3784cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3794cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3804cbe5319SToby Isaac                 v[d] = 0.;
381ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3824cbe5319SToby Isaac               }
383ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = v[d];
3844cbe5319SToby Isaac             }
3854cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3864cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3874cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3884cbe5319SToby Isaac                   v[d] = 0.;
389ad540459SPierre Jolivet                   for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3904cbe5319SToby Isaac                 }
391ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = v[d];
3924cbe5319SToby Isaac               }
3934cbe5319SToby Isaac             }
3944cbe5319SToby Isaac           }
395a6216909SToby Isaac         }
396a6216909SToby Isaac       }
3976858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
3989566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&T));
399a6216909SToby Isaac     }
400a6216909SToby Isaac     offset += nq;
401a6216909SToby Isaac   }
402a6216909SToby Isaac   {
403a6216909SToby Isaac     MPI_Datatype origtype;
4044cbe5319SToby Isaac 
405a6216909SToby Isaac     if (datatype == PETSC_SCALAR) {
406a6216909SToby Isaac       origtype = MPIU_SCALAR;
407a6216909SToby Isaac     } else {
408a6216909SToby Isaac       origtype = MPIU_REAL;
409a6216909SToby Isaac     }
410a6216909SToby Isaac     if (B) {
411a6216909SToby Isaac       MPI_Datatype Btype;
412a6216909SToby Isaac 
4136497c311SBarry Smith       PetscCallMPI(MPI_Type_contiguous((PetscMPIInt)nc, origtype, &Btype));
4149566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Btype));
4159566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4169566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4179566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Btype));
418a6216909SToby Isaac     }
419a6216909SToby Isaac     if (D) {
420a6216909SToby Isaac       MPI_Datatype Dtype;
421a6216909SToby Isaac 
4226497c311SBarry Smith       PetscCallMPI(MPI_Type_contiguous((PetscMPIInt)(nc * dim), origtype, &Dtype));
4239566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Dtype));
4249566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4259566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4269566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Dtype));
427a6216909SToby Isaac     }
428a6216909SToby Isaac     if (H) {
429a6216909SToby Isaac       MPI_Datatype Htype;
430a6216909SToby Isaac 
4316497c311SBarry Smith       PetscCallMPI(MPI_Type_contiguous((PetscMPIInt)(nc * dim * dim), origtype, &Htype));
4329566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Htype));
4339566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4349566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4359566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Htype));
436a6216909SToby Isaac     }
437a6216909SToby Isaac   }
4389566063dSJacob Faibussowitsch   PetscCall(PetscFree4(v, J, invJ, detJ));
4399566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBr, cellDr, cellHr));
4409566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBs, cellDs, cellHs));
4419566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellPoints, coordsReal, coordsRef));
4429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&pointType));
4439566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
445a6216909SToby Isaac }
446a6216909SToby Isaac 
447d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H)
448d71ae5a4SJacob Faibussowitsch {
449805e7170SToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
450805e7170SToby Isaac   PetscInt         h, imin;
451805e7170SToby Isaac   PetscInt         dim;
452805e7170SToby Isaac   PetscClassId     id;
453805e7170SToby Isaac   PetscQuadrature  quad = NULL;
454b7260050SToby Isaac   PetscInt         maxDegree;
455805e7170SToby Isaac   PetscFEGeom     *geom;
456805e7170SToby Isaac   PetscInt         Nq, Nc, dimC, qNc, N;
457805e7170SToby Isaac   PetscInt         numPoints;
458805e7170SToby Isaac   void            *qB = NULL, *qD = NULL, *qH = NULL;
459805e7170SToby Isaac   const PetscReal *weights;
460805e7170SToby Isaac   MPI_Datatype     mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL;
461805e7170SToby Isaac   PetscObject      disc;
462805e7170SToby Isaac   DMField          coordField;
463805e7170SToby Isaac 
464805e7170SToby Isaac   PetscFunctionBegin;
465805e7170SToby Isaac   Nc = field->numComponents;
4669566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dimC));
4679566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dim));
4689566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numPoints));
4699566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, NULL));
470805e7170SToby Isaac   for (h = 0; h < dsfield->height; h++) {
471805e7170SToby Isaac     PetscInt hEnd;
472805e7170SToby Isaac 
4739566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
474805e7170SToby Isaac     if (imin < hEnd) break;
475805e7170SToby Isaac   }
476805e7170SToby Isaac   dim -= h;
4776858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
4789566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
47908401ef6SPierre Jolivet   PetscCheck(id == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported");
4809566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(field->dm, &coordField));
4819566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
48248a46eb9SPierre Jolivet   if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
4839566063dSJacob Faibussowitsch   if (!quad) PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, &quad));
48428b400f6SJacob Faibussowitsch   PetscCheck(quad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages");
4859566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom));
4869566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights));
48708401ef6SPierre Jolivet   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components");
488805e7170SToby Isaac   N = numPoints * Nq * Nc;
4899566063dSJacob Faibussowitsch   if (B) PetscCall(DMGetWorkArray(field->dm, N, mpitype, &qB));
4909566063dSJacob Faibussowitsch   if (D) PetscCall(DMGetWorkArray(field->dm, N * dimC, mpitype, &qD));
4919566063dSJacob Faibussowitsch   if (H) PetscCall(DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
4929566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, pointIS, quad, type, qB, qD, qH));
493805e7170SToby Isaac   if (B) {
494805e7170SToby Isaac     PetscInt i, j, k;
495805e7170SToby Isaac 
496805e7170SToby Isaac     if (type == PETSC_SCALAR) {
497805e7170SToby Isaac       PetscScalar *sB  = (PetscScalar *)B;
498805e7170SToby Isaac       PetscScalar *sqB = (PetscScalar *)qB;
499805e7170SToby Isaac 
500805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
501805e7170SToby Isaac         PetscReal vol = 0.;
502805e7170SToby Isaac 
503ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) sB[i * Nc + j] = 0.;
504805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
505805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
506ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[(i * Nq + k) * Nc + j];
507805e7170SToby Isaac         }
50862bd480fSMatthew G. Knepley         for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol;
509805e7170SToby Isaac       }
510805e7170SToby Isaac     } else {
511805e7170SToby Isaac       PetscReal *rB  = (PetscReal *)B;
512805e7170SToby Isaac       PetscReal *rqB = (PetscReal *)qB;
513805e7170SToby Isaac 
514805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
515805e7170SToby Isaac         PetscReal vol = 0.;
516805e7170SToby Isaac 
517ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) rB[i * Nc + j] = 0.;
518805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
519805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
520ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) rB[i * Nc + j] += weights[k] * rqB[(i * Nq + k) * Nc + j];
521805e7170SToby Isaac         }
522805e7170SToby Isaac         for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol;
523805e7170SToby Isaac       }
524805e7170SToby Isaac     }
525805e7170SToby Isaac   }
526805e7170SToby Isaac   if (D) {
527805e7170SToby Isaac     PetscInt i, j, k, l, m;
528805e7170SToby Isaac 
529805e7170SToby Isaac     if (type == PETSC_SCALAR) {
530805e7170SToby Isaac       PetscScalar *sD  = (PetscScalar *)D;
531805e7170SToby Isaac       PetscScalar *sqD = (PetscScalar *)qD;
532805e7170SToby Isaac 
533805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
534805e7170SToby Isaac         PetscReal vol = 0.;
535805e7170SToby Isaac 
536ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) sD[i * Nc * dimC + j] = 0.;
537805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
538805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
539805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
540805e7170SToby Isaac             PetscScalar pD[3] = {0., 0., 0.};
541805e7170SToby Isaac 
542805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
543ad540459SPierre 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];
544805e7170SToby Isaac             }
545ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
546805e7170SToby Isaac           }
547805e7170SToby Isaac         }
548805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol;
549805e7170SToby Isaac       }
550805e7170SToby Isaac     } else {
551805e7170SToby Isaac       PetscReal *rD  = (PetscReal *)D;
552805e7170SToby Isaac       PetscReal *rqD = (PetscReal *)qD;
553805e7170SToby Isaac 
554805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
555805e7170SToby Isaac         PetscReal vol = 0.;
556805e7170SToby Isaac 
557ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) rD[i * Nc * dimC + j] = 0.;
558805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
559805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
560805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
561805e7170SToby Isaac             PetscReal pD[3] = {0., 0., 0.};
562805e7170SToby Isaac 
563805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
564ad540459SPierre 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];
565805e7170SToby Isaac             }
566ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
567805e7170SToby Isaac           }
568805e7170SToby Isaac         }
569805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol;
570805e7170SToby Isaac       }
571805e7170SToby Isaac     }
572805e7170SToby Isaac   }
573805e7170SToby Isaac   if (H) {
574805e7170SToby Isaac     PetscInt i, j, k, l, m, q, r;
575805e7170SToby Isaac 
576805e7170SToby Isaac     if (type == PETSC_SCALAR) {
577805e7170SToby Isaac       PetscScalar *sH  = (PetscScalar *)H;
578805e7170SToby Isaac       PetscScalar *sqH = (PetscScalar *)qH;
579805e7170SToby Isaac 
580805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
581805e7170SToby Isaac         PetscReal vol = 0.;
582805e7170SToby Isaac 
583ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) sH[i * Nc * dimC * dimC + j] = 0.;
584805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
5854d1a973fSToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
586805e7170SToby Isaac 
587805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
588805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
5899371c9d4SSatish Balay             PetscScalar pH[3][3] = {
5909371c9d4SSatish Balay               {0., 0., 0.},
5919371c9d4SSatish Balay               {0., 0., 0.},
5929371c9d4SSatish Balay               {0., 0., 0.}
5939371c9d4SSatish Balay             };
594805e7170SToby Isaac             const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC];
595805e7170SToby Isaac 
596805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
597805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
598805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
599ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r];
600805e7170SToby Isaac                 }
601805e7170SToby Isaac               }
602805e7170SToby Isaac             }
603805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
604ad540459SPierre 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];
605805e7170SToby Isaac             }
606805e7170SToby Isaac           }
607805e7170SToby Isaac         }
608805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol;
609805e7170SToby Isaac       }
610805e7170SToby Isaac     } else {
611805e7170SToby Isaac       PetscReal *rH  = (PetscReal *)H;
612805e7170SToby Isaac       PetscReal *rqH = (PetscReal *)qH;
613805e7170SToby Isaac 
614805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
615805e7170SToby Isaac         PetscReal vol = 0.;
616805e7170SToby Isaac 
617ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) rH[i * Nc * dimC * dimC + j] = 0.;
618805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
619805e7170SToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
620805e7170SToby Isaac 
621805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
622805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
6239371c9d4SSatish Balay             PetscReal pH[3][3] = {
6249371c9d4SSatish Balay               {0., 0., 0.},
6259371c9d4SSatish Balay               {0., 0., 0.},
6269371c9d4SSatish Balay               {0., 0., 0.}
6279371c9d4SSatish Balay             };
628805e7170SToby Isaac             const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC];
629805e7170SToby Isaac 
630805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
631805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
632805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
633ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r];
634805e7170SToby Isaac                 }
635805e7170SToby Isaac               }
636805e7170SToby Isaac             }
637805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
638ad540459SPierre 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];
639805e7170SToby Isaac             }
640805e7170SToby Isaac           }
641805e7170SToby Isaac         }
642805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol;
643805e7170SToby Isaac       }
644805e7170SToby Isaac     }
645805e7170SToby Isaac   }
6469566063dSJacob Faibussowitsch   if (B) PetscCall(DMRestoreWorkArray(field->dm, N, mpitype, &qB));
6479566063dSJacob Faibussowitsch   if (D) PetscCall(DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD));
6489566063dSJacob Faibussowitsch   if (H) PetscCall(DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
6499566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
6509566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
6513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
652805e7170SToby Isaac }
653805e7170SToby Isaac 
654d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
655d71ae5a4SJacob Faibussowitsch {
65651a454edSToby Isaac   DMField_DS  *dsfield;
65751a454edSToby Isaac   PetscObject  disc;
658741db25dSToby Isaac   PetscInt     h, imin, imax;
65951a454edSToby Isaac   PetscClassId id;
66051a454edSToby Isaac 
66151a454edSToby Isaac   PetscFunctionBegin;
66251a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
6639566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
664741db25dSToby Isaac   if (imin >= imax) {
66551aa0bf7SToby Isaac     h = 0;
66651aa0bf7SToby Isaac   } else {
667f99c8401SToby Isaac     for (h = 0; h < dsfield->height; h++) {
66851a454edSToby Isaac       PetscInt hEnd;
66951a454edSToby Isaac 
6709566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
67151a454edSToby Isaac       if (imin < hEnd) break;
67251a454edSToby Isaac     }
67351aa0bf7SToby Isaac   }
6746858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
67651a454edSToby Isaac   if (id == PETSCFE_CLASSID) {
67751a454edSToby Isaac     PetscFE    fe = (PetscFE)disc;
67851a454edSToby Isaac     PetscSpace sp;
67951a454edSToby Isaac 
6809566063dSJacob Faibussowitsch     PetscCall(PetscFEGetBasisSpace(fe, &sp));
6819566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetDegree(sp, minDegree, maxDegree));
68251a454edSToby Isaac   }
6833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68451a454edSToby Isaac }
68551a454edSToby Isaac 
686d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetFVQuadrature_Internal(DMField field, IS pointIS, PetscQuadrature *quad)
687d71ae5a4SJacob Faibussowitsch {
688a2aba86cSMatthew G. Knepley   DM              dm = field->dm;
689a2aba86cSMatthew G. Knepley   const PetscInt *points;
690a2aba86cSMatthew G. Knepley   DMPolytopeType  ct;
691a2aba86cSMatthew G. Knepley   PetscInt        dim, n;
692a2aba86cSMatthew G. Knepley   PetscBool       isplex;
693a2aba86cSMatthew G. Knepley 
694a2aba86cSMatthew G. Knepley   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex));
6969566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &n));
697a2aba86cSMatthew G. Knepley   if (isplex && n) {
6989566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
6999566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pointIS, &points));
7009566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, points[0], &ct));
701a2aba86cSMatthew G. Knepley     switch (ct) {
702a2aba86cSMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
703d71ae5a4SJacob Faibussowitsch     case DM_POLYTOPE_TETRAHEDRON:
704d71ae5a4SJacob Faibussowitsch       PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 1, -1.0, 1.0, quad));
705d71ae5a4SJacob Faibussowitsch       break;
706d71ae5a4SJacob Faibussowitsch     default:
707d71ae5a4SJacob Faibussowitsch       PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
708a2aba86cSMatthew G. Knepley     }
7099566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pointIS, &points));
7101baa6e33SBarry Smith   } else PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, quad));
7113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
712a2aba86cSMatthew G. Knepley }
713a2aba86cSMatthew G. Knepley 
714d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad)
715d71ae5a4SJacob Faibussowitsch {
7161fdf1f07SMatthew G. Knepley   PetscInt     h, dim, imax, imin, cellHeight;
71751a454edSToby Isaac   DM           dm;
71851a454edSToby Isaac   DMField_DS  *dsfield;
71951a454edSToby Isaac   PetscObject  disc;
72051a454edSToby Isaac   PetscFE      fe;
72151a454edSToby Isaac   PetscClassId id;
72251a454edSToby Isaac 
72351a454edSToby Isaac   PetscFunctionBegin;
72451a454edSToby Isaac   dm      = field->dm;
72551a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
7269566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
7279566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
72851a454edSToby Isaac   for (h = 0; h <= dim; h++) {
72951a454edSToby Isaac     PetscInt hStart, hEnd;
73051a454edSToby Isaac 
7319566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
7326c041b70SSatish Balay     if (imax >= hStart && imin < hEnd) break;
73351a454edSToby Isaac   }
7349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
7351fdf1f07SMatthew G. Knepley   h -= cellHeight;
73651a454edSToby Isaac   *quad = NULL;
737f99c8401SToby Isaac   if (h < dsfield->height) {
7386858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
7399566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
7403ba16761SJacob Faibussowitsch     if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
74151a454edSToby Isaac     fe = (PetscFE)disc;
7429566063dSJacob Faibussowitsch     PetscCall(PetscFEGetQuadrature(fe, quad));
7439566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)*quad));
74451a454edSToby Isaac   }
7453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74651a454edSToby Isaac }
74751a454edSToby Isaac 
748d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom)
749d71ae5a4SJacob Faibussowitsch {
7500145028aSToby Isaac   const PetscInt *points;
7510145028aSToby Isaac   PetscInt        p, dim, dE, numFaces, Nq;
752b7260050SToby Isaac   PetscInt        maxDegree;
7530145028aSToby Isaac   DMLabel         depthLabel;
7540145028aSToby Isaac   IS              cellIS;
7550145028aSToby Isaac   DM              dm = field->dm;
7560145028aSToby Isaac 
7570145028aSToby Isaac   PetscFunctionBegin;
7580145028aSToby Isaac   dim = geom->dim;
7590145028aSToby Isaac   dE  = geom->dimEmbed;
7609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
7619566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS));
7629566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(field, cellIS, NULL, &maxDegree));
7639566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointIS, &points));
7640145028aSToby Isaac   numFaces = geom->numCells;
7650145028aSToby Isaac   Nq       = geom->numPoints;
766f15274beSMatthew Knepley   /* First, set local faces and flip normals so that they are outward for the first supporting cell */
767f15274beSMatthew Knepley   for (p = 0; p < numFaces; p++) {
768f15274beSMatthew Knepley     PetscInt        point = points[p];
769f15274beSMatthew Knepley     PetscInt        suppSize, s, coneSize, c, numChildren;
7700e18dc48SMatthew G. Knepley     const PetscInt *supp;
771f15274beSMatthew Knepley 
7729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
77328b400f6SJacob Faibussowitsch     PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
7749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, point, &suppSize));
77563a3b9bcSJacob Faibussowitsch     PetscCheck(suppSize <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, suppSize);
776f15274beSMatthew Knepley     if (!suppSize) continue;
7779566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, point, &supp));
778f15274beSMatthew Knepley     for (s = 0; s < suppSize; ++s) {
7790e18dc48SMatthew G. Knepley       const PetscInt *cone, *ornt;
7800e18dc48SMatthew G. Knepley 
7819566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, supp[s], &coneSize));
7820e18dc48SMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, supp[s], &cone, &ornt));
7839371c9d4SSatish Balay       for (c = 0; c < coneSize; ++c)
7849371c9d4SSatish Balay         if (cone[c] == point) break;
7850e18dc48SMatthew G. Knepley       geom->face[p][s * 2 + 0] = c;
7860e18dc48SMatthew G. Knepley       geom->face[p][s * 2 + 1] = ornt[c];
78784cbd072SMatthew G. Knepley       PetscCall(DMPlexRestoreOrientedCone(dm, supp[s], &cone, &ornt));
78884cbd072SMatthew G. Knepley       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]);
789f15274beSMatthew Knepley     }
7900e18dc48SMatthew G. Knepley     if (geom->face[p][1] < 0) {
791f15274beSMatthew Knepley       PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d;
792f15274beSMatthew Knepley 
7939371c9d4SSatish Balay       for (q = 0; q < Np; ++q)
7949371c9d4SSatish Balay         for (d = 0; d < dE; ++d) geom->n[(p * Np + q) * dE + d] = -geom->n[(p * Np + q) * dE + d];
795f15274beSMatthew Knepley     }
796f15274beSMatthew Knepley   }
797b7260050SToby Isaac   if (maxDegree <= 1) {
798313f16f2SMatthew G. Knepley     PetscQuadrature cellQuad = NULL;
7990145028aSToby Isaac     PetscInt        numCells, offset, *cells;
8000145028aSToby Isaac     PetscFEGeom    *cellGeom;
8010145028aSToby Isaac     IS              suppIS;
8020145028aSToby Isaac 
803313f16f2SMatthew G. Knepley     if (quad) {
804313f16f2SMatthew G. Knepley       DM         dm;
805313f16f2SMatthew G. Knepley       PetscReal *points, *weights;
806313f16f2SMatthew G. Knepley       PetscInt   tdim, Nc, Np;
807313f16f2SMatthew G. Knepley 
808313f16f2SMatthew G. Knepley       PetscCall(DMFieldGetDM(field, &dm));
809313f16f2SMatthew G. Knepley       PetscCall(DMGetDimension(dm, &tdim));
810313f16f2SMatthew G. Knepley       if (tdim > dim) {
811313f16f2SMatthew G. Knepley         // Make a compatible cell quadrature (points don't matter since its affine)
812313f16f2SMatthew G. Knepley         PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad));
813313f16f2SMatthew G. Knepley         PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Np, NULL, NULL));
814313f16f2SMatthew G. Knepley         PetscCall(PetscCalloc1((dim + 1) * Np, &points));
815313f16f2SMatthew G. Knepley         PetscCall(PetscCalloc1(Nc * Np, &weights));
816313f16f2SMatthew G. Knepley         PetscCall(PetscQuadratureSetData(cellQuad, dim + 1, Nc, Np, points, weights));
817313f16f2SMatthew G. Knepley       } else {
818313f16f2SMatthew G. Knepley         // TODO J will be wrong here, but other things need to be fixed
819313f16f2SMatthew G. Knepley         //   This path comes from calling DMProjectBdFieldLabelLocal() in Plex ex5
820313f16f2SMatthew G. Knepley         PetscCall(PetscObjectReference((PetscObject)quad));
821313f16f2SMatthew G. Knepley         cellQuad = quad;
822313f16f2SMatthew G. Knepley       }
823313f16f2SMatthew G. Knepley     }
8240145028aSToby Isaac     for (p = 0, numCells = 0; p < numFaces; p++) {
8250145028aSToby Isaac       PetscInt point = points[p];
8260145028aSToby Isaac       PetscInt numSupp, numChildren;
8270145028aSToby Isaac 
8289566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
82928b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
8309566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
83163a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
8320145028aSToby Isaac       numCells += numSupp;
8330145028aSToby Isaac     }
8349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numCells, &cells));
8350145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8360145028aSToby Isaac       PetscInt        point = points[p];
8370145028aSToby Isaac       PetscInt        numSupp, s;
8380145028aSToby Isaac       const PetscInt *supp;
8390145028aSToby Isaac 
8409566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8419566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
842ad540459SPierre Jolivet       for (s = 0; s < numSupp; s++, offset++) cells[offset] = supp[s];
8430145028aSToby Isaac     }
8449566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
845313f16f2SMatthew G. Knepley     PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FALSE, &cellGeom));
8460145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8470145028aSToby Isaac       PetscInt        point = points[p];
8480145028aSToby Isaac       PetscInt        numSupp, s, q;
8490145028aSToby Isaac       const PetscInt *supp;
8500145028aSToby Isaac 
8519566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8529566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
8530145028aSToby Isaac       for (s = 0; s < numSupp; s++, offset++) {
8540145028aSToby Isaac         for (q = 0; q < Nq * dE * dE; q++) {
85541418987SMatthew G. Knepley           geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
8560145028aSToby Isaac           geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
8570145028aSToby Isaac         }
85841418987SMatthew G. Knepley         for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
8590145028aSToby Isaac       }
8600145028aSToby Isaac     }
8619566063dSJacob Faibussowitsch     PetscCall(PetscFEGeomDestroy(&cellGeom));
862313f16f2SMatthew G. Knepley     PetscCall(PetscQuadratureDestroy(&cellQuad));
8639566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&suppIS));
8649566063dSJacob Faibussowitsch     PetscCall(PetscFree(cells));
8650145028aSToby Isaac   } else {
8666858538eSMatthew G. Knepley     DMField_DS    *dsfield = (DMField_DS *)field->data;
8670145028aSToby Isaac     PetscObject    faceDisc, cellDisc;
8680145028aSToby Isaac     PetscClassId   faceId, cellId;
8690145028aSToby Isaac     PetscDualSpace dsp;
8700145028aSToby Isaac     DM             K;
871ba2698f1SMatthew G. Knepley     DMPolytopeType ct;
8720145028aSToby Isaac     PetscInt(*co)[2][3];
8730145028aSToby Isaac     PetscInt        coneSize;
8740145028aSToby Isaac     PetscInt      **counts;
8750145028aSToby Isaac     PetscInt        f, i, o, q, s;
876e59eb9b3SMatthew G. Knepley     PetscBool       found = PETSC_FALSE;
8770145028aSToby Isaac     const PetscInt *coneK;
878ba2698f1SMatthew G. Knepley     PetscInt        eStart, minOrient, maxOrient, numOrient;
8790145028aSToby Isaac     PetscInt       *orients;
8800145028aSToby Isaac     PetscReal     **orientPoints;
8810145028aSToby Isaac     PetscReal      *cellPoints;
8820145028aSToby Isaac     PetscReal      *dummyWeights;
8830145028aSToby Isaac     PetscQuadrature cellQuad = NULL;
8840145028aSToby Isaac 
8856858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 1, dsfield->disc, &faceDisc));
8866858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
8879566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(faceDisc, &faceId));
8889566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(cellDisc, &cellId));
8891dca8a05SBarry Smith     PetscCheck(faceId == PETSCFE_CLASSID && cellId == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported");
8909566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace((PetscFE)cellDisc, &dsp));
8919566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(dsp, &K));
8929566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(K, 1, &eStart, NULL));
8939566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(K, eStart, &ct));
8949566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(K, 0, &coneSize));
8959566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
8969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numFaces, &co, coneSize, &counts));
8979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dE * Nq, &cellPoints));
8989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nq, &dummyWeights));
8999566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad));
9009566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights));
9010145028aSToby Isaac     minOrient = PETSC_MAX_INT;
9020145028aSToby Isaac     maxOrient = PETSC_MIN_INT;
9030145028aSToby Isaac     for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */
9040145028aSToby Isaac       PetscInt        point = points[p];
9050145028aSToby Isaac       PetscInt        numSupp, numChildren;
9060145028aSToby Isaac       const PetscInt *supp;
9070145028aSToby Isaac 
9089566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
90928b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
9109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
91163a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
9129566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
9130145028aSToby Isaac       for (s = 0; s < numSupp; s++) {
9140145028aSToby Isaac         PetscInt        cell = supp[s];
9150145028aSToby Isaac         PetscInt        numCone;
9160145028aSToby Isaac         const PetscInt *cone, *orient;
9170145028aSToby Isaac 
9189566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(dm, cell, &numCone));
919e59eb9b3SMatthew G. Knepley         // When we extract submeshes, we hang cells from the side that are not fully realized. We ignore these
920e59eb9b3SMatthew G. Knepley         if (numCone == 1) {
921e59eb9b3SMatthew G. Knepley           co[p][s][0] = -1;
922e59eb9b3SMatthew G. Knepley           co[p][s][1] = -1;
923e59eb9b3SMatthew G. Knepley           co[p][s][2] = -1;
924e59eb9b3SMatthew G. Knepley           continue;
925e59eb9b3SMatthew G. Knepley         }
92608401ef6SPierre Jolivet         PetscCheck(numCone == coneSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element");
9279566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(dm, cell, &cone));
9289566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeOrientation(dm, cell, &orient));
9290145028aSToby Isaac         for (f = 0; f < coneSize; f++) {
9300145028aSToby Isaac           if (cone[f] == point) break;
9310145028aSToby Isaac         }
9320145028aSToby Isaac         co[p][s][0] = f;
9330145028aSToby Isaac         co[p][s][1] = orient[f];
9340145028aSToby Isaac         co[p][s][2] = cell;
9350145028aSToby Isaac         minOrient   = PetscMin(minOrient, orient[f]);
936cd93b0e1SMatthew G. Knepley         maxOrient   = PetscMax(maxOrient, orient[f]);
937e59eb9b3SMatthew G. Knepley         found       = PETSC_TRUE;
9380145028aSToby Isaac       }
9390145028aSToby Isaac       for (; s < 2; s++) {
9400145028aSToby Isaac         co[p][s][0] = -1;
9410145028aSToby Isaac         co[p][s][1] = -1;
9420145028aSToby Isaac         co[p][s][2] = -1;
9430145028aSToby Isaac       }
9440145028aSToby Isaac     }
945e59eb9b3SMatthew G. Knepley     numOrient = found ? maxOrient + 1 - minOrient : 0;
9469566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
9470145028aSToby Isaac     /* count all (face,orientation) doubles that appear */
9489566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(numOrient, &orients, numOrient, &orientPoints));
9499566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscCalloc1(numOrient + 1, &counts[f]));
9500145028aSToby Isaac     for (p = 0; p < numFaces; p++) {
9510145028aSToby Isaac       for (s = 0; s < 2; s++) {
9520145028aSToby Isaac         if (co[p][s][0] >= 0) {
9530145028aSToby Isaac           counts[co[p][s][0]][co[p][s][1] - minOrient]++;
9540145028aSToby Isaac           orients[co[p][s][1] - minOrient]++;
9550145028aSToby Isaac         }
9560145028aSToby Isaac       }
9570145028aSToby Isaac     }
9580145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
9590145028aSToby Isaac       if (orients[o]) {
9600145028aSToby Isaac         PetscInt orient = o + minOrient;
9610145028aSToby Isaac         PetscInt q;
9620145028aSToby Isaac 
9639566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Nq * dim, &orientPoints[o]));
9640145028aSToby Isaac         /* rotate the quadrature points appropriately */
965ba2698f1SMatthew G. Knepley         switch (ct) {
966d71ae5a4SJacob Faibussowitsch         case DM_POLYTOPE_POINT:
967d71ae5a4SJacob Faibussowitsch           break;
968ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_SEGMENT:
9690145028aSToby Isaac           if (orient == -2 || orient == 1) {
970ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = -geom->xi[q];
9710145028aSToby Isaac           } else {
972ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = geom->xi[q];
9730145028aSToby Isaac           }
9740145028aSToby Isaac           break;
975ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:
9760145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9770145028aSToby Isaac             PetscReal lambda[3];
9780145028aSToby Isaac             PetscReal lambdao[3];
9790145028aSToby Isaac 
9800145028aSToby Isaac             /* convert to barycentric */
9810145028aSToby Isaac             lambda[0] = -(geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.;
9820145028aSToby Isaac             lambda[1] = (geom->xi[2 * q] + 1.) / 2.;
9830145028aSToby Isaac             lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.;
9840145028aSToby Isaac             if (orient >= 0) {
985ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(orient + i) % 3];
9860145028aSToby Isaac             } else {
987ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(-(orient + i) + 3) % 3];
9880145028aSToby Isaac             }
9890145028aSToby Isaac             /* convert to coordinates */
9900145028aSToby Isaac             orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1];
9910145028aSToby Isaac             orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2];
9920145028aSToby Isaac           }
9930145028aSToby Isaac           break;
994ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL:
9950145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9960145028aSToby Isaac             PetscReal xi[2], xio[2];
9970145028aSToby Isaac             PetscInt  oabs = (orient >= 0) ? orient : -(orient + 1);
9980145028aSToby Isaac 
9990145028aSToby Isaac             xi[0] = geom->xi[2 * q];
10000145028aSToby Isaac             xi[1] = geom->xi[2 * q + 1];
10010145028aSToby Isaac             switch (oabs) {
10020145028aSToby Isaac             case 1:
10030145028aSToby Isaac               xio[0] = xi[1];
10040145028aSToby Isaac               xio[1] = -xi[0];
10050145028aSToby Isaac               break;
1006f4d061e9SPierre Jolivet             case 2:
1007f4d061e9SPierre Jolivet               xio[0] = -xi[0];
1008f4d061e9SPierre Jolivet               xio[1] = -xi[1];
1009f4d061e9SPierre Jolivet               break;
1010f4d061e9SPierre Jolivet             case 3:
1011f4d061e9SPierre Jolivet               xio[0] = -xi[1];
1012f4d061e9SPierre Jolivet               xio[1] = xi[0];
1013f4d061e9SPierre Jolivet               break;
1014aa1371cbSToby Isaac             case 0:
1015aa1371cbSToby Isaac             default:
1016aa1371cbSToby Isaac               xio[0] = xi[0];
1017aa1371cbSToby Isaac               xio[1] = xi[1];
1018aa1371cbSToby Isaac               break;
10190145028aSToby Isaac             }
1020ad540459SPierre Jolivet             if (orient < 0) xio[0] = -xio[0];
10210145028aSToby Isaac             orientPoints[o][2 * q + 0] = xio[0];
10220145028aSToby Isaac             orientPoints[o][2 * q + 1] = xio[1];
10230145028aSToby Isaac           }
10240145028aSToby Isaac           break;
1025d71ae5a4SJacob Faibussowitsch         default:
1026d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s not yet supported", DMPolytopeTypes[ct]);
10270145028aSToby Isaac         }
10280145028aSToby Isaac       }
10290145028aSToby Isaac     }
10300145028aSToby Isaac     for (f = 0; f < coneSize; f++) {
10310145028aSToby Isaac       PetscInt  face = coneK[f];
10324d1a973fSToby Isaac       PetscReal v0[3];
10332d89661fSMatthew G. Knepley       PetscReal J[9], detJ;
10340145028aSToby Isaac       PetscInt  numCells, offset;
10350145028aSToby Isaac       PetscInt *cells;
10360145028aSToby Isaac       IS        suppIS;
10370145028aSToby Isaac 
10389566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ));
10390145028aSToby Isaac       for (o = 0; o <= numOrient; o++) {
10400145028aSToby Isaac         PetscFEGeom *cellGeom;
10410145028aSToby Isaac 
10420145028aSToby Isaac         if (!counts[f][o]) continue;
10430145028aSToby Isaac         /* If this (face,orientation) double appears,
10440145028aSToby Isaac          * convert the face quadrature points into volume quadrature points */
10450145028aSToby Isaac         for (q = 0; q < Nq; q++) {
10460145028aSToby Isaac           PetscReal xi0[3] = {-1., -1., -1.};
10470145028aSToby Isaac 
10480145028aSToby Isaac           CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]);
10490145028aSToby Isaac         }
10500145028aSToby Isaac         for (p = 0, numCells = 0; p < numFaces; p++) {
10510145028aSToby Isaac           for (s = 0; s < 2; s++) {
10520145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++;
10530145028aSToby Isaac           }
10540145028aSToby Isaac         }
10559566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(numCells, &cells));
10560145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10570145028aSToby Isaac           for (s = 0; s < 2; s++) {
1058ad540459SPierre Jolivet             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) cells[offset++] = co[p][s][2];
10590145028aSToby Isaac           }
10600145028aSToby Isaac         }
10619566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
10629566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FALSE, &cellGeom));
10630145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10640145028aSToby Isaac           for (s = 0; s < 2; s++) {
10650145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
10660145028aSToby Isaac               for (q = 0; q < Nq * dE * dE; q++) {
106741418987SMatthew G. Knepley                 geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
10680145028aSToby Isaac                 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
10690145028aSToby Isaac               }
107041418987SMatthew G. Knepley               for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
10710145028aSToby Isaac               offset++;
10720145028aSToby Isaac             }
10730145028aSToby Isaac           }
10740145028aSToby Isaac         }
10759566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&cellGeom));
10769566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&suppIS));
10779566063dSJacob Faibussowitsch         PetscCall(PetscFree(cells));
10780145028aSToby Isaac       }
10790145028aSToby Isaac     }
10800145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
108148a46eb9SPierre Jolivet       if (orients[o]) PetscCall(PetscFree(orientPoints[o]));
10820145028aSToby Isaac     }
10839566063dSJacob Faibussowitsch     PetscCall(PetscFree2(orients, orientPoints));
10849566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&cellQuad));
10859566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscFree(counts[f]));
10869566063dSJacob Faibussowitsch     PetscCall(PetscFree2(co, counts));
10870145028aSToby Isaac   }
10889566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointIS, &points));
10899566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
10903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10910145028aSToby Isaac }
10920145028aSToby Isaac 
1093d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldInitialize_DS(DMField field)
1094d71ae5a4SJacob Faibussowitsch {
109551a454edSToby Isaac   PetscFunctionBegin;
109651a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DS;
109751a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DS;
109851a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DS;
1099805e7170SToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DS;
1100b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DS;
110151a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS;
110251a454edSToby Isaac   field->ops->view                    = DMFieldView_DS;
11030145028aSToby Isaac   field->ops->computeFaceData         = DMFieldComputeFaceData_DS;
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110551a454edSToby Isaac }
110651a454edSToby Isaac 
1107d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field)
1108d71ae5a4SJacob Faibussowitsch {
110951a454edSToby Isaac   DMField_DS *dsfield;
111051a454edSToby Isaac 
111151a454edSToby Isaac   PetscFunctionBegin;
11124dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dsfield));
111351a454edSToby Isaac   field->data = dsfield;
11149566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_DS(field));
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111651a454edSToby Isaac }
111751a454edSToby Isaac 
1118d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDSWithDG(DM dm, DM dmDG, PetscInt fieldNum, Vec vec, Vec vecDG, DMField *field)
1119d71ae5a4SJacob Faibussowitsch {
112051a454edSToby Isaac   DMField      b;
112151a454edSToby Isaac   DMField_DS  *dsfield;
11226858538eSMatthew G. Knepley   PetscObject  disc = NULL, discDG = NULL;
11236858538eSMatthew G. Knepley   PetscSection section;
112451a454edSToby Isaac   PetscBool    isContainer   = PETSC_FALSE;
112551a454edSToby Isaac   PetscClassId id            = -1;
112651a454edSToby Isaac   PetscInt     numComponents = -1, dsNumFields;
112751a454edSToby Isaac 
112851a454edSToby Isaac   PetscFunctionBegin;
11299566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
11309566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetFieldComponents(section, fieldNum, &numComponents));
11319566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &dsNumFields));
11329566063dSJacob Faibussowitsch   if (dsNumFields) PetscCall(DMGetField(dm, fieldNum, NULL, &disc));
11336858538eSMatthew G. Knepley   if (dsNumFields && dmDG) {
11346858538eSMatthew G. Knepley     PetscCall(DMGetField(dmDG, fieldNum, NULL, &discDG));
11356858538eSMatthew G. Knepley     PetscCall(PetscObjectReference(discDG));
11366858538eSMatthew G. Knepley   }
113751a454edSToby Isaac   if (disc) {
11389566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
113951a454edSToby Isaac     isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE;
114051a454edSToby Isaac   }
114151a454edSToby Isaac   if (!disc || isContainer) {
114251a454edSToby Isaac     MPI_Comm       comm = PetscObjectComm((PetscObject)dm);
114351a454edSToby Isaac     PetscFE        fe;
11442df84da0SMatthew G. Knepley     DMPolytopeType ct, locct = DM_POLYTOPE_UNKNOWN;
11452df84da0SMatthew G. Knepley     PetscInt       dim, cStart, cEnd, cellHeight;
114651a454edSToby Isaac 
11479566063dSJacob Faibussowitsch     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
11489566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
11499566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
11509566063dSJacob Faibussowitsch     if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &locct));
1151*6a210b70SBarry Smith     PetscCallMPI(MPIU_Allreduce(&locct, &ct, 1, MPI_INT, MPI_MIN, comm));
11529566063dSJacob Faibussowitsch     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, numComponents, ct, 1, PETSC_DETERMINE, &fe));
11539566063dSJacob Faibussowitsch     PetscCall(PetscFEViewFromOptions(fe, NULL, "-field_fe_view"));
115451a454edSToby Isaac     disc = (PetscObject)fe;
11551baa6e33SBarry Smith   } else PetscCall(PetscObjectReference(disc));
11569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
11576858538eSMatthew G. Knepley   if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)disc, &numComponents));
11586858538eSMatthew G. Knepley   else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine number of discretization components");
11599566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm, numComponents, DMFIELD_VERTEX, &b));
11609566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b, DMFIELDDS));
116151a454edSToby Isaac   dsfield           = (DMField_DS *)b->data;
116251a454edSToby Isaac   dsfield->fieldNum = fieldNum;
11639566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dsfield->height));
1164f99c8401SToby Isaac   dsfield->height++;
11659566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(dsfield->height, &dsfield->disc));
116651a454edSToby Isaac   dsfield->disc[0] = disc;
11679566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)vec));
116851a454edSToby Isaac   dsfield->vec = vec;
11696858538eSMatthew G. Knepley   if (dmDG) {
11706858538eSMatthew G. Knepley     dsfield->dmDG = dmDG;
11716858538eSMatthew G. Knepley     PetscCall(PetscCalloc1(dsfield->height, &dsfield->discDG));
11726858538eSMatthew G. Knepley     dsfield->discDG[0] = discDG;
11736858538eSMatthew G. Knepley     PetscCall(PetscObjectReference((PetscObject)vecDG));
11746858538eSMatthew G. Knepley     dsfield->vecDG = vecDG;
11756858538eSMatthew G. Knepley   }
117651a454edSToby Isaac   *field = b;
11773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11786858538eSMatthew G. Knepley }
117951a454edSToby Isaac 
1180d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec, DMField *field)
1181d71ae5a4SJacob Faibussowitsch {
11826858538eSMatthew G. Knepley   PetscFunctionBegin;
11836858538eSMatthew G. Knepley   PetscCall(DMFieldCreateDSWithDG(dm, NULL, fieldNum, vec, NULL, field));
11843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118551a454edSToby Isaac }
1186