xref: /petsc/src/dm/field/impls/ds/dmfieldds.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
189371c9d4SSatish Balay static PetscErrorCode DMFieldDestroy_DS(DMField field) {
196858538eSMatthew G. Knepley   DMField_DS *dsfield = (DMField_DS *)field->data;
2051a454edSToby Isaac   PetscInt    i;
2151a454edSToby Isaac 
2251a454edSToby Isaac   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dsfield->vec));
246858538eSMatthew G. Knepley   for (i = 0; i < dsfield->height; i++) PetscCall(PetscObjectDereference(dsfield->disc[i]));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(dsfield->disc));
266858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dsfield->vecDG));
279371c9d4SSatish Balay   if (dsfield->discDG)
289371c9d4SSatish Balay     for (i = 0; i < dsfield->height; i++) PetscCall(PetscObjectDereference(dsfield->discDG[i]));
296858538eSMatthew G. Knepley   PetscCall(PetscFree(dsfield->discDG));
309566063dSJacob Faibussowitsch   PetscCall(PetscFree(dsfield));
3151a454edSToby Isaac   PetscFunctionReturn(0);
3251a454edSToby Isaac }
3351a454edSToby Isaac 
349371c9d4SSatish Balay static PetscErrorCode DMFieldView_DS(DMField field, PetscViewer viewer) {
3551a454edSToby Isaac   DMField_DS *dsfield = (DMField_DS *)field->data;
3651a454edSToby Isaac   PetscObject disc;
376858538eSMatthew G. Knepley   PetscBool   iascii;
3851a454edSToby Isaac 
3951a454edSToby Isaac   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4151a454edSToby Isaac   disc = dsfield->disc[0];
4251a454edSToby Isaac   if (iascii) {
4363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PetscDS field %" PetscInt_FMT "\n", dsfield->fieldNum));
449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
459566063dSJacob Faibussowitsch     PetscCall(PetscObjectView(disc, viewer));
469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
4751a454edSToby Isaac   }
489566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
4928b400f6SJacob Faibussowitsch   PetscCheck(!dsfield->multifieldVec, PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "View of subfield not implemented yet");
509566063dSJacob Faibussowitsch   PetscCall(VecView(dsfield->vec, viewer));
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
5251a454edSToby Isaac   PetscFunctionReturn(0);
5351a454edSToby Isaac }
5451a454edSToby Isaac 
559371c9d4SSatish Balay static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject discList[], PetscObject *disc) {
5651a454edSToby Isaac   PetscFunctionBegin;
576858538eSMatthew G. Knepley   if (!discList[height]) {
5851a454edSToby Isaac     PetscClassId id;
5951a454edSToby Isaac 
606858538eSMatthew G. Knepley     PetscCall(PetscObjectGetClassId(discList[0], &id));
616858538eSMatthew G. Knepley     if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateHeightTrace((PetscFE)discList[0], height, (PetscFE *)&discList[height]));
6251a454edSToby Isaac   }
636858538eSMatthew G. Knepley   *disc = discList[height];
6451a454edSToby Isaac   PetscFunctionReturn(0);
6551a454edSToby Isaac }
6651a454edSToby Isaac 
6762bd480fSMatthew G. Knepley /* y[m,c] = A[m,n,c] . b[n] */
6851a454edSToby Isaac #define DMFieldDSdot(y, A, b, m, n, c, cast) \
6951a454edSToby Isaac   do { \
7051a454edSToby Isaac     PetscInt _i, _j, _k; \
7151a454edSToby Isaac     for (_i = 0; _i < (m); _i++) { \
72ad540459SPierre Jolivet       for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] = 0.; \
7351a454edSToby Isaac       for (_j = 0; _j < (n); _j++) { \
74ad540459SPierre Jolivet         for (_k = 0; _k < (c); _k++) (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \
7551a454edSToby Isaac       } \
7651a454edSToby Isaac     } \
7751a454edSToby Isaac   } while (0)
7851a454edSToby Isaac 
796858538eSMatthew G. Knepley /*
806858538eSMatthew 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().
816858538eSMatthew G. Knepley */
829371c9d4SSatish Balay PetscErrorCode DMFieldGetClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[]) {
836858538eSMatthew G. Knepley   DMField_DS        *dsfield = (DMField_DS *)field->data;
846858538eSMatthew G. Knepley   DM                 fdm     = dsfield->dmDG;
856858538eSMatthew G. Knepley   PetscSection       s       = NULL;
866858538eSMatthew G. Knepley   const PetscScalar *cvalues;
876858538eSMatthew G. Knepley   PetscInt           pStart, pEnd;
886858538eSMatthew G. Knepley 
896858538eSMatthew G. Knepley   PetscFunctionBeginHot;
906858538eSMatthew G. Knepley   *isDG   = PETSC_FALSE;
916858538eSMatthew G. Knepley   *Nc     = 0;
926858538eSMatthew G. Knepley   *array  = NULL;
936858538eSMatthew G. Knepley   *values = NULL;
946858538eSMatthew G. Knepley   /* Check for cellwise section */
956858538eSMatthew G. Knepley   if (fdm) PetscCall(DMGetLocalSection(fdm, &s));
966858538eSMatthew G. Knepley   if (!s) goto cg;
976858538eSMatthew G. Knepley   /* Check that the cell exists in the cellwise section */
986858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
996858538eSMatthew G. Knepley   if (cell < pStart || cell >= pEnd) goto cg;
1006858538eSMatthew G. Knepley   /* Check for cellwise coordinates for this cell */
1016858538eSMatthew G. Knepley   PetscCall(PetscSectionGetDof(s, cell, Nc));
1026858538eSMatthew G. Knepley   if (!*Nc) goto cg;
1036858538eSMatthew G. Knepley   /* Check for cellwise coordinates */
1046858538eSMatthew G. Knepley   if (!dsfield->vecDG) goto cg;
1056858538eSMatthew G. Knepley   /* Get cellwise coordinates */
1066858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(dsfield->vecDG, array));
1076858538eSMatthew G. Knepley   PetscCall(DMPlexPointLocalRead(fdm, cell, *array, &cvalues));
1086858538eSMatthew G. Knepley   PetscCall(DMGetWorkArray(fdm, *Nc, MPIU_SCALAR, values));
1096858538eSMatthew G. Knepley   PetscCall(PetscArraycpy(*values, cvalues, *Nc));
1106858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(dsfield->vecDG, array));
1116858538eSMatthew G. Knepley   *isDG = PETSC_TRUE;
1126858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1136858538eSMatthew G. Knepley cg:
1146858538eSMatthew G. Knepley   /* Use continuous values */
1156858538eSMatthew G. Knepley   PetscCall(DMFieldGetDM(field, &fdm));
1166858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(fdm, &s));
1176267c18cSMatthew G. Knepley   PetscCall(PetscSectionGetField(s, dsfield->fieldNum, &s));
1186858538eSMatthew G. Knepley   PetscCall(DMPlexVecGetClosure(fdm, s, dsfield->vec, cell, Nc, values));
1196858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1206858538eSMatthew G. Knepley }
1216858538eSMatthew G. Knepley 
1229371c9d4SSatish Balay PetscErrorCode DMFieldRestoreClosure_Internal(DMField field, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *values[]) {
1236858538eSMatthew G. Knepley   DMField_DS  *dsfield = (DMField_DS *)field->data;
1246858538eSMatthew G. Knepley   DM           fdm;
1256858538eSMatthew G. Knepley   PetscSection s;
1266858538eSMatthew G. Knepley 
1276858538eSMatthew G. Knepley   PetscFunctionBeginHot;
1286858538eSMatthew G. Knepley   if (*isDG) {
1296858538eSMatthew G. Knepley     PetscCall(DMRestoreWorkArray(dsfield->dmDG, *Nc, MPIU_SCALAR, values));
1306858538eSMatthew G. Knepley   } else {
1316858538eSMatthew G. Knepley     PetscCall(DMFieldGetDM(field, &fdm));
1326858538eSMatthew G. Knepley     PetscCall(DMGetLocalSection(fdm, &s));
1336858538eSMatthew G. Knepley     PetscCall(DMPlexVecRestoreClosure(fdm, s, dsfield->vec, cell, Nc, (PetscScalar **)values));
1346858538eSMatthew G. Knepley   }
1356858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1366858538eSMatthew G. Knepley }
1376858538eSMatthew G. Knepley 
138ef0bb6c7SMatthew G. Knepley /* TODO: Reorganize interface so that I can reuse a tabulation rather than mallocing each time */
1399371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) {
14051a454edSToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
14151a454edSToby Isaac   DM               dm;
14251a454edSToby Isaac   PetscObject      disc;
14351a454edSToby Isaac   PetscClassId     classid;
14451a454edSToby Isaac   PetscInt         nq, nc, dim, meshDim, numCells;
14551a454edSToby Isaac   PetscSection     section;
14651a454edSToby Isaac   const PetscReal *qpoints;
14751a454edSToby Isaac   PetscBool        isStride;
14851a454edSToby Isaac   const PetscInt  *points = NULL;
14951a454edSToby Isaac   PetscInt         sfirst = -1, stride = -1;
15051a454edSToby Isaac 
15151a454edSToby Isaac   PetscFunctionBeginHot;
15251a454edSToby Isaac   dm = field->dm;
15351a454edSToby Isaac   nc = field->numComponents;
1549566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &nq, &qpoints, NULL));
1556858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, dsfield->height - 1 - dim, dsfield->disc, &disc));
1569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &meshDim));
1579566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1589566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetField(section, dsfield->fieldNum, &section));
1599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &classid));
16051a454edSToby Isaac   /* TODO: batch */
1619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
1629566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
1631baa6e33SBarry Smith   if (isStride) PetscCall(ISStrideGetInfo(pointIS, &sfirst, &stride));
1641baa6e33SBarry Smith   else PetscCall(ISGetIndices(pointIS, &points));
16551a454edSToby Isaac   if (classid == PETSCFE_CLASSID) {
16651a454edSToby Isaac     PetscFE         fe = (PetscFE)disc;
16751a454edSToby Isaac     PetscInt        feDim, i;
168ef0bb6c7SMatthew G. Knepley     PetscInt        K = H ? 2 : (D ? 1 : (B ? 0 : -1));
169ef0bb6c7SMatthew G. Knepley     PetscTabulation T;
17051a454edSToby Isaac 
1719566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDimension(fe, &feDim));
1729566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fe, 1, nq, qpoints, K, &T));
17351a454edSToby Isaac     for (i = 0; i < numCells; i++) {
17451a454edSToby Isaac       PetscInt           c = isStride ? (sfirst + i * stride) : points[i];
1752710589cSToby Isaac       PetscInt           closureSize;
1766858538eSMatthew G. Knepley       const PetscScalar *array;
1772710589cSToby Isaac       PetscScalar       *elem = NULL;
1786858538eSMatthew G. Knepley       PetscBool          isDG;
17951a454edSToby Isaac 
1806858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
18151a454edSToby Isaac       if (B) {
18262bd480fSMatthew G. Knepley         /* field[c] = T[q,b,c] . coef[b], so v[c] = T[q,b,c] . coords[b] */
18351a454edSToby Isaac         if (type == PETSC_SCALAR) {
18451a454edSToby Isaac           PetscScalar *cB = &((PetscScalar *)B)[nc * nq * i];
18551a454edSToby Isaac 
186ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
18751a454edSToby Isaac         } else {
18851a454edSToby Isaac           PetscReal *cB = &((PetscReal *)B)[nc * nq * i];
18951a454edSToby Isaac 
190ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
19151a454edSToby Isaac         }
19251a454edSToby Isaac       }
19351a454edSToby Isaac       if (D) {
19451a454edSToby Isaac         if (type == PETSC_SCALAR) {
19551a454edSToby Isaac           PetscScalar *cD = &((PetscScalar *)D)[nc * nq * dim * i];
19651a454edSToby Isaac 
197ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
19851a454edSToby Isaac         } else {
19951a454edSToby Isaac           PetscReal *cD = &((PetscReal *)D)[nc * nq * dim * i];
20051a454edSToby Isaac 
201ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
20251a454edSToby Isaac         }
20351a454edSToby Isaac       }
20451a454edSToby Isaac       if (H) {
20551a454edSToby Isaac         if (type == PETSC_SCALAR) {
20651a454edSToby Isaac           PetscScalar *cH = &((PetscScalar *)H)[nc * nq * dim * dim * i];
20751a454edSToby Isaac 
208ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
20951a454edSToby Isaac         } else {
21051a454edSToby Isaac           PetscReal *cH = &((PetscReal *)H)[nc * nq * dim * dim * i];
21151a454edSToby Isaac 
212ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
21351a454edSToby Isaac         }
21451a454edSToby Isaac       }
2156858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
21651a454edSToby Isaac     }
2179566063dSJacob Faibussowitsch     PetscCall(PetscTabulationDestroy(&T));
2182da392ccSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented");
21948a46eb9SPierre Jolivet   if (!isStride) PetscCall(ISRestoreIndices(pointIS, &points));
22051a454edSToby Isaac   PetscFunctionReturn(0);
22151a454edSToby Isaac }
22251a454edSToby Isaac 
2239371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) {
224a6216909SToby Isaac   DMField_DS        *dsfield = (DMField_DS *)field->data;
225a6216909SToby Isaac   PetscSF            cellSF  = NULL;
226a6216909SToby Isaac   const PetscSFNode *cells;
227a6216909SToby Isaac   PetscInt           c, nFound, numCells, feDim, nc;
228a6216909SToby Isaac   const PetscInt    *cellDegrees;
229a6216909SToby Isaac   const PetscScalar *pointsArray;
230a6216909SToby Isaac   PetscScalar       *cellPoints;
231a6216909SToby Isaac   PetscInt           gatherSize, gatherMax;
232a6216909SToby Isaac   PetscInt           dim, dimR, offset;
233a6216909SToby Isaac   MPI_Datatype       pointType;
234a6216909SToby Isaac   PetscObject        cellDisc;
235a6216909SToby Isaac   PetscFE            cellFE;
236a6216909SToby Isaac   PetscClassId       discID;
237a6216909SToby Isaac   PetscReal         *coordsReal, *coordsRef;
238a6216909SToby Isaac   PetscSection       section;
239a6216909SToby Isaac   PetscScalar       *cellBs = NULL, *cellDs = NULL, *cellHs = NULL;
240a6216909SToby Isaac   PetscReal         *cellBr = NULL, *cellDr = NULL, *cellHr = NULL;
2414cbe5319SToby Isaac   PetscReal         *v, *J, *invJ, *detJ;
242a6216909SToby Isaac 
243a6216909SToby Isaac   PetscFunctionBegin;
244a6216909SToby Isaac   nc = field->numComponents;
2459566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(field->dm, &section));
2466858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(cellDisc, &discID));
24808401ef6SPierre Jolivet   PetscCheck(discID == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization type not supported");
249a6216909SToby Isaac   cellFE = (PetscFE)cellDisc;
2509566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDimension(cellFE, &feDim));
2519566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dim));
2529566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dimR));
2539566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF));
2549566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells));
255ad540459SPierre 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);
2569566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(cellSF, &cellDegrees));
2579566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(cellSF, &cellDegrees));
258a6216909SToby Isaac   for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) {
259a6216909SToby Isaac     gatherMax = PetscMax(gatherMax, cellDegrees[c]);
260a6216909SToby Isaac     gatherSize += cellDegrees[c];
261a6216909SToby Isaac   }
2629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(gatherSize * dim, &cellPoints, gatherMax * dim, &coordsReal, gatherMax * dimR, &coordsRef));
2639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(gatherMax * dimR, &v, gatherMax * dimR * dimR, &J, gatherMax * dimR * dimR, &invJ, gatherMax, &detJ));
2641baa6e33SBarry 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));
2651baa6e33SBarry Smith   else PetscCall(PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr));
266a6216909SToby Isaac 
2679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(dim, MPIU_SCALAR, &pointType));
2689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&pointType));
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(points, &pointsArray));
2709566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints));
2719566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints));
2729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(points, &pointsArray));
273a6216909SToby Isaac   for (c = 0, offset = 0; c < numCells; c++) {
274a6216909SToby Isaac     PetscInt nq = cellDegrees[c], p;
275a6216909SToby Isaac 
276a6216909SToby Isaac     if (nq) {
277ef0bb6c7SMatthew G. Knepley       PetscInt           K = H ? 2 : (D ? 1 : (B ? 0 : -1));
278ef0bb6c7SMatthew G. Knepley       PetscTabulation    T;
2796858538eSMatthew G. Knepley       PetscQuadrature    quad;
2806858538eSMatthew G. Knepley       const PetscScalar *array;
281a6216909SToby Isaac       PetscScalar       *elem = NULL;
2824cbe5319SToby Isaac       PetscReal         *quadPoints;
2836858538eSMatthew G. Knepley       PetscBool          isDG;
2846858538eSMatthew G. Knepley       PetscInt           closureSize, d, e, f, g;
285a6216909SToby Isaac 
2864d1a973fSToby Isaac       for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]);
2879566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef));
2889566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(cellFE, 1, nq, coordsRef, K, &T));
2899566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad));
2909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(dimR * nq, &quadPoints));
2914cbe5319SToby Isaac       for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p];
2929566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL));
2939566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ));
2949566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
2956858538eSMatthew G. Knepley       PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
296a6216909SToby Isaac       if (B) {
297a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
298a6216909SToby Isaac           PetscScalar *cB = &cellBs[nc * offset];
299a6216909SToby Isaac 
300ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar));
301a6216909SToby Isaac         } else {
302a6216909SToby Isaac           PetscReal *cB = &cellBr[nc * offset];
303a6216909SToby Isaac 
304ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart);
305a6216909SToby Isaac         }
306a6216909SToby Isaac       }
307a6216909SToby Isaac       if (D) {
308a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
309a6216909SToby Isaac           PetscScalar *cD = &cellDs[nc * dim * offset];
310a6216909SToby Isaac 
311ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar));
3124cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3134cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3144d1a973fSToby Isaac               PetscScalar vs[3];
3154d1a973fSToby Isaac 
3164cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3174d1a973fSToby Isaac                 vs[d] = 0.;
318ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3194cbe5319SToby Isaac               }
320ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = vs[d];
3214cbe5319SToby Isaac             }
3224cbe5319SToby Isaac           }
323a6216909SToby Isaac         } else {
324a6216909SToby Isaac           PetscReal *cD = &cellDr[nc * dim * offset];
325a6216909SToby Isaac 
326ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart);
3274cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3284cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3294cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3304cbe5319SToby Isaac                 v[d] = 0.;
331ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
3324cbe5319SToby Isaac               }
333ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = v[d];
3344cbe5319SToby Isaac             }
3354cbe5319SToby Isaac           }
336a6216909SToby Isaac         }
337a6216909SToby Isaac       }
338a6216909SToby Isaac       if (H) {
339a6216909SToby Isaac         if (datatype == PETSC_SCALAR) {
340a6216909SToby Isaac           PetscScalar *cH = &cellHs[nc * dim * dim * offset];
341a6216909SToby Isaac 
342ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar));
3434cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3444cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3454d1a973fSToby Isaac               PetscScalar vs[3];
3464d1a973fSToby Isaac 
3474cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3484d1a973fSToby Isaac                 vs[d] = 0.;
349ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3504cbe5319SToby Isaac               }
351ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = vs[d];
3524cbe5319SToby Isaac             }
3534cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3544cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
3554d1a973fSToby Isaac                 PetscScalar vs[3];
3564d1a973fSToby Isaac 
3574cbe5319SToby Isaac                 for (d = 0; d < dimR; d++) {
3584d1a973fSToby Isaac                   vs[d] = 0.;
359ad540459SPierre Jolivet                   for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
3604cbe5319SToby Isaac                 }
361ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d];
3624cbe5319SToby Isaac               }
3634cbe5319SToby Isaac             }
3644cbe5319SToby Isaac           }
365a6216909SToby Isaac         } else {
366a6216909SToby Isaac           PetscReal *cH = &cellHr[nc * dim * dim * offset];
367a6216909SToby Isaac 
368ef0bb6c7SMatthew G. Knepley           DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart);
3694cbe5319SToby Isaac           for (p = 0; p < nq; p++) {
3704cbe5319SToby Isaac             for (g = 0; g < nc * dimR; g++) {
3714cbe5319SToby Isaac               for (d = 0; d < dimR; d++) {
3724cbe5319SToby Isaac                 v[d] = 0.;
373ad540459SPierre Jolivet                 for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
3744cbe5319SToby Isaac               }
375ad540459SPierre Jolivet               for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = v[d];
3764cbe5319SToby Isaac             }
3774cbe5319SToby Isaac             for (g = 0; g < nc; g++) {
3784cbe5319SToby Isaac               for (f = 0; f < dimR; f++) {
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 * p + g) * dimR + e) * dimR + f];
3824cbe5319SToby Isaac                 }
383ad540459SPierre Jolivet                 for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = v[d];
3844cbe5319SToby Isaac               }
3854cbe5319SToby Isaac             }
3864cbe5319SToby Isaac           }
387a6216909SToby Isaac         }
388a6216909SToby Isaac       }
3896858538eSMatthew G. Knepley       PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem));
3909566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&T));
391a6216909SToby Isaac     }
392a6216909SToby Isaac     offset += nq;
393a6216909SToby Isaac   }
394a6216909SToby Isaac   {
395a6216909SToby Isaac     MPI_Datatype origtype;
3964cbe5319SToby Isaac 
397a6216909SToby Isaac     if (datatype == PETSC_SCALAR) {
398a6216909SToby Isaac       origtype = MPIU_SCALAR;
399a6216909SToby Isaac     } else {
400a6216909SToby Isaac       origtype = MPIU_REAL;
401a6216909SToby Isaac     }
402a6216909SToby Isaac     if (B) {
403a6216909SToby Isaac       MPI_Datatype Btype;
404a6216909SToby Isaac 
4059566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc, origtype, &Btype));
4069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Btype));
4079566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4089566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B));
4099566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Btype));
410a6216909SToby Isaac     }
411a6216909SToby Isaac     if (D) {
412a6216909SToby Isaac       MPI_Datatype Dtype;
413a6216909SToby Isaac 
4149566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc * dim, origtype, &Dtype));
4159566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Dtype));
4169566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4179566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D));
4189566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Dtype));
419a6216909SToby Isaac     }
420a6216909SToby Isaac     if (H) {
421a6216909SToby Isaac       MPI_Datatype Htype;
422a6216909SToby Isaac 
4239566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_contiguous(nc * dim * dim, origtype, &Htype));
4249566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_commit(&Htype));
4259566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterBegin(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4269566063dSJacob Faibussowitsch       PetscCall(PetscSFScatterEnd(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H));
4279566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Type_free(&Htype));
428a6216909SToby Isaac     }
429a6216909SToby Isaac   }
4309566063dSJacob Faibussowitsch   PetscCall(PetscFree4(v, J, invJ, detJ));
4319566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBr, cellDr, cellHr));
4329566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellBs, cellDs, cellHs));
4339566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellPoints, coordsReal, coordsRef));
4349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&pointType));
4359566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
436a6216909SToby Isaac   PetscFunctionReturn(0);
437a6216909SToby Isaac }
438a6216909SToby Isaac 
4399371c9d4SSatish Balay static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) {
440805e7170SToby Isaac   DMField_DS      *dsfield = (DMField_DS *)field->data;
441805e7170SToby Isaac   PetscInt         h, imin;
442805e7170SToby Isaac   PetscInt         dim;
443805e7170SToby Isaac   PetscClassId     id;
444805e7170SToby Isaac   PetscQuadrature  quad = NULL;
445b7260050SToby Isaac   PetscInt         maxDegree;
446805e7170SToby Isaac   PetscFEGeom     *geom;
447805e7170SToby Isaac   PetscInt         Nq, Nc, dimC, qNc, N;
448805e7170SToby Isaac   PetscInt         numPoints;
449805e7170SToby Isaac   void            *qB = NULL, *qD = NULL, *qH = NULL;
450805e7170SToby Isaac   const PetscReal *weights;
451805e7170SToby Isaac   MPI_Datatype     mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL;
452805e7170SToby Isaac   PetscObject      disc;
453805e7170SToby Isaac   DMField          coordField;
454805e7170SToby Isaac 
455805e7170SToby Isaac   PetscFunctionBegin;
456805e7170SToby Isaac   Nc = field->numComponents;
4579566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(field->dm, &dimC));
4589566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(field->dm, &dim));
4599566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numPoints));
4609566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, NULL));
461805e7170SToby Isaac   for (h = 0; h < dsfield->height; h++) {
462805e7170SToby Isaac     PetscInt hEnd;
463805e7170SToby Isaac 
4649566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
465805e7170SToby Isaac     if (imin < hEnd) break;
466805e7170SToby Isaac   }
467805e7170SToby Isaac   dim -= h;
4686858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
47008401ef6SPierre Jolivet   PetscCheck(id == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported");
4719566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(field->dm, &coordField));
4729566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
47348a46eb9SPierre Jolivet   if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
4749566063dSJacob Faibussowitsch   if (!quad) PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, &quad));
47528b400f6SJacob Faibussowitsch   PetscCheck(quad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages");
4769566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom));
4779566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights));
47808401ef6SPierre Jolivet   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components");
479805e7170SToby Isaac   N = numPoints * Nq * Nc;
4809566063dSJacob Faibussowitsch   if (B) PetscCall(DMGetWorkArray(field->dm, N, mpitype, &qB));
4819566063dSJacob Faibussowitsch   if (D) PetscCall(DMGetWorkArray(field->dm, N * dimC, mpitype, &qD));
4829566063dSJacob Faibussowitsch   if (H) PetscCall(DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
4839566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, pointIS, quad, type, qB, qD, qH));
484805e7170SToby Isaac   if (B) {
485805e7170SToby Isaac     PetscInt i, j, k;
486805e7170SToby Isaac 
487805e7170SToby Isaac     if (type == PETSC_SCALAR) {
488805e7170SToby Isaac       PetscScalar *sB  = (PetscScalar *)B;
489805e7170SToby Isaac       PetscScalar *sqB = (PetscScalar *)qB;
490805e7170SToby Isaac 
491805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
492805e7170SToby Isaac         PetscReal vol = 0.;
493805e7170SToby Isaac 
494ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) sB[i * Nc + j] = 0.;
495805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
496805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
497ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[(i * Nq + k) * Nc + j];
498805e7170SToby Isaac         }
49962bd480fSMatthew G. Knepley         for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol;
500805e7170SToby Isaac       }
501805e7170SToby Isaac     } else {
502805e7170SToby Isaac       PetscReal *rB  = (PetscReal *)B;
503805e7170SToby Isaac       PetscReal *rqB = (PetscReal *)qB;
504805e7170SToby Isaac 
505805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
506805e7170SToby Isaac         PetscReal vol = 0.;
507805e7170SToby Isaac 
508ad540459SPierre Jolivet         for (j = 0; j < Nc; j++) rB[i * Nc + j] = 0.;
509805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
510805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
511ad540459SPierre Jolivet           for (j = 0; j < Nc; j++) rB[i * Nc + j] += weights[k] * rqB[(i * Nq + k) * Nc + j];
512805e7170SToby Isaac         }
513805e7170SToby Isaac         for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol;
514805e7170SToby Isaac       }
515805e7170SToby Isaac     }
516805e7170SToby Isaac   }
517805e7170SToby Isaac   if (D) {
518805e7170SToby Isaac     PetscInt i, j, k, l, m;
519805e7170SToby Isaac 
520805e7170SToby Isaac     if (type == PETSC_SCALAR) {
521805e7170SToby Isaac       PetscScalar *sD  = (PetscScalar *)D;
522805e7170SToby Isaac       PetscScalar *sqD = (PetscScalar *)qD;
523805e7170SToby Isaac 
524805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
525805e7170SToby Isaac         PetscReal vol = 0.;
526805e7170SToby Isaac 
527ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) sD[i * Nc * dimC + j] = 0.;
528805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
529805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
530805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
531805e7170SToby Isaac             PetscScalar pD[3] = {0., 0., 0.};
532805e7170SToby Isaac 
533805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
534ad540459SPierre 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];
535805e7170SToby Isaac             }
536ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
537805e7170SToby Isaac           }
538805e7170SToby Isaac         }
539805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol;
540805e7170SToby Isaac       }
541805e7170SToby Isaac     } else {
542805e7170SToby Isaac       PetscReal *rD  = (PetscReal *)D;
543805e7170SToby Isaac       PetscReal *rqD = (PetscReal *)qD;
544805e7170SToby Isaac 
545805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
546805e7170SToby Isaac         PetscReal vol = 0.;
547805e7170SToby Isaac 
548ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC; j++) rD[i * Nc * dimC + j] = 0.;
549805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
550805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
551805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
552805e7170SToby Isaac             PetscReal pD[3] = {0., 0., 0.};
553805e7170SToby Isaac 
554805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
555ad540459SPierre 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];
556805e7170SToby Isaac             }
557ad540459SPierre Jolivet             for (l = 0; l < dimC; l++) rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
558805e7170SToby Isaac           }
559805e7170SToby Isaac         }
560805e7170SToby Isaac         for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol;
561805e7170SToby Isaac       }
562805e7170SToby Isaac     }
563805e7170SToby Isaac   }
564805e7170SToby Isaac   if (H) {
565805e7170SToby Isaac     PetscInt i, j, k, l, m, q, r;
566805e7170SToby Isaac 
567805e7170SToby Isaac     if (type == PETSC_SCALAR) {
568805e7170SToby Isaac       PetscScalar *sH  = (PetscScalar *)H;
569805e7170SToby Isaac       PetscScalar *sqH = (PetscScalar *)qH;
570805e7170SToby Isaac 
571805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
572805e7170SToby Isaac         PetscReal vol = 0.;
573805e7170SToby Isaac 
574ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) sH[i * Nc * dimC * dimC + j] = 0.;
575805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
5764d1a973fSToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
577805e7170SToby Isaac 
578805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
579805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
5809371c9d4SSatish Balay             PetscScalar pH[3][3] = {
5819371c9d4SSatish Balay               {0., 0., 0.},
5829371c9d4SSatish Balay               {0., 0., 0.},
5839371c9d4SSatish Balay               {0., 0., 0.}
5849371c9d4SSatish Balay             };
585805e7170SToby Isaac             const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC];
586805e7170SToby Isaac 
587805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
588805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
589805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
590ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r];
591805e7170SToby Isaac                 }
592805e7170SToby Isaac               }
593805e7170SToby Isaac             }
594805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
595ad540459SPierre 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];
596805e7170SToby Isaac             }
597805e7170SToby Isaac           }
598805e7170SToby Isaac         }
599805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol;
600805e7170SToby Isaac       }
601805e7170SToby Isaac     } else {
602805e7170SToby Isaac       PetscReal *rH  = (PetscReal *)H;
603805e7170SToby Isaac       PetscReal *rqH = (PetscReal *)qH;
604805e7170SToby Isaac 
605805e7170SToby Isaac       for (i = 0; i < numPoints; i++) {
606805e7170SToby Isaac         PetscReal vol = 0.;
607805e7170SToby Isaac 
608ad540459SPierre Jolivet         for (j = 0; j < Nc * dimC * dimC; j++) rH[i * Nc * dimC * dimC + j] = 0.;
609805e7170SToby Isaac         for (k = 0; k < Nq; k++) {
610805e7170SToby Isaac           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
611805e7170SToby Isaac 
612805e7170SToby Isaac           vol += geom->detJ[i * Nq + k] * weights[k];
613805e7170SToby Isaac           for (j = 0; j < Nc; j++) {
6149371c9d4SSatish Balay             PetscReal pH[3][3] = {
6159371c9d4SSatish Balay               {0., 0., 0.},
6169371c9d4SSatish Balay               {0., 0., 0.},
6179371c9d4SSatish Balay               {0., 0., 0.}
6189371c9d4SSatish Balay             };
619805e7170SToby Isaac             const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC];
620805e7170SToby Isaac 
621805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
622805e7170SToby Isaac               for (m = 0; m < dimC; m++) {
623805e7170SToby Isaac                 for (q = 0; q < dim; q++) {
624ad540459SPierre Jolivet                   for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r];
625805e7170SToby Isaac                 }
626805e7170SToby Isaac               }
627805e7170SToby Isaac             }
628805e7170SToby Isaac             for (l = 0; l < dimC; l++) {
629ad540459SPierre 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];
630805e7170SToby Isaac             }
631805e7170SToby Isaac           }
632805e7170SToby Isaac         }
633805e7170SToby Isaac         for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol;
634805e7170SToby Isaac       }
635805e7170SToby Isaac     }
636805e7170SToby Isaac   }
6379566063dSJacob Faibussowitsch   if (B) PetscCall(DMRestoreWorkArray(field->dm, N, mpitype, &qB));
6389566063dSJacob Faibussowitsch   if (D) PetscCall(DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD));
6399566063dSJacob Faibussowitsch   if (H) PetscCall(DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH));
6409566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
6419566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
642805e7170SToby Isaac   PetscFunctionReturn(0);
643805e7170SToby Isaac }
644805e7170SToby Isaac 
6459371c9d4SSatish Balay static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree) {
64651a454edSToby Isaac   DMField_DS  *dsfield;
64751a454edSToby Isaac   PetscObject  disc;
648741db25dSToby Isaac   PetscInt     h, imin, imax;
64951a454edSToby Isaac   PetscClassId id;
65051a454edSToby Isaac 
65151a454edSToby Isaac   PetscFunctionBegin;
65251a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
6539566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
654741db25dSToby Isaac   if (imin >= imax) {
65551aa0bf7SToby Isaac     h = 0;
65651aa0bf7SToby Isaac   } else {
657f99c8401SToby Isaac     for (h = 0; h < dsfield->height; h++) {
65851a454edSToby Isaac       PetscInt hEnd;
65951a454edSToby Isaac 
6609566063dSJacob Faibussowitsch       PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd));
66151a454edSToby Isaac       if (imin < hEnd) break;
66251a454edSToby Isaac     }
66351aa0bf7SToby Isaac   }
6646858538eSMatthew G. Knepley   PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
6659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
66651a454edSToby Isaac   if (id == PETSCFE_CLASSID) {
66751a454edSToby Isaac     PetscFE    fe = (PetscFE)disc;
66851a454edSToby Isaac     PetscSpace sp;
66951a454edSToby Isaac 
6709566063dSJacob Faibussowitsch     PetscCall(PetscFEGetBasisSpace(fe, &sp));
6719566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetDegree(sp, minDegree, maxDegree));
67251a454edSToby Isaac   }
67351a454edSToby Isaac   PetscFunctionReturn(0);
67451a454edSToby Isaac }
67551a454edSToby Isaac 
6769371c9d4SSatish Balay PetscErrorCode DMFieldGetFVQuadrature_Internal(DMField field, IS pointIS, PetscQuadrature *quad) {
677a2aba86cSMatthew G. Knepley   DM              dm = field->dm;
678a2aba86cSMatthew G. Knepley   const PetscInt *points;
679a2aba86cSMatthew G. Knepley   DMPolytopeType  ct;
680a2aba86cSMatthew G. Knepley   PetscInt        dim, n;
681a2aba86cSMatthew G. Knepley   PetscBool       isplex;
682a2aba86cSMatthew G. Knepley 
683a2aba86cSMatthew G. Knepley   PetscFunctionBegin;
6849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex));
6859566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &n));
686a2aba86cSMatthew G. Knepley   if (isplex && n) {
6879566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
6889566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pointIS, &points));
6899566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, points[0], &ct));
690a2aba86cSMatthew G. Knepley     switch (ct) {
691a2aba86cSMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
6929371c9d4SSatish Balay     case DM_POLYTOPE_TETRAHEDRON: PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 1, -1.0, 1.0, quad)); break;
6939566063dSJacob Faibussowitsch     default: PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
694a2aba86cSMatthew G. Knepley     }
6959566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pointIS, &points));
6961baa6e33SBarry Smith   } else PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, quad));
697a2aba86cSMatthew G. Knepley   PetscFunctionReturn(0);
698a2aba86cSMatthew G. Knepley }
699a2aba86cSMatthew G. Knepley 
7009371c9d4SSatish Balay static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) {
7011fdf1f07SMatthew G. Knepley   PetscInt     h, dim, imax, imin, cellHeight;
70251a454edSToby Isaac   DM           dm;
70351a454edSToby Isaac   DMField_DS  *dsfield;
70451a454edSToby Isaac   PetscObject  disc;
70551a454edSToby Isaac   PetscFE      fe;
70651a454edSToby Isaac   PetscClassId id;
70751a454edSToby Isaac 
70851a454edSToby Isaac   PetscFunctionBegin;
70951a454edSToby Isaac   dm      = field->dm;
71051a454edSToby Isaac   dsfield = (DMField_DS *)field->data;
7119566063dSJacob Faibussowitsch   PetscCall(ISGetMinMax(pointIS, &imin, &imax));
7129566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
71351a454edSToby Isaac   for (h = 0; h <= dim; h++) {
71451a454edSToby Isaac     PetscInt hStart, hEnd;
71551a454edSToby Isaac 
7169566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
7176c041b70SSatish Balay     if (imax >= hStart && imin < hEnd) break;
71851a454edSToby Isaac   }
7199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
7201fdf1f07SMatthew G. Knepley   h -= cellHeight;
72151a454edSToby Isaac   *quad = NULL;
722f99c8401SToby Isaac   if (h < dsfield->height) {
7236858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc));
7249566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
72551a454edSToby Isaac     if (id != PETSCFE_CLASSID) PetscFunctionReturn(0);
72651a454edSToby Isaac     fe = (PetscFE)disc;
7279566063dSJacob Faibussowitsch     PetscCall(PetscFEGetQuadrature(fe, quad));
7289566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)*quad));
72951a454edSToby Isaac   }
73051a454edSToby Isaac   PetscFunctionReturn(0);
73151a454edSToby Isaac }
73251a454edSToby Isaac 
7339371c9d4SSatish Balay static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom) {
7340145028aSToby Isaac   const PetscInt *points;
7350145028aSToby Isaac   PetscInt        p, dim, dE, numFaces, Nq;
736b7260050SToby Isaac   PetscInt        maxDegree;
7370145028aSToby Isaac   DMLabel         depthLabel;
7380145028aSToby Isaac   IS              cellIS;
7390145028aSToby Isaac   DM              dm = field->dm;
7400145028aSToby Isaac 
7410145028aSToby Isaac   PetscFunctionBegin;
7420145028aSToby Isaac   dim = geom->dim;
7430145028aSToby Isaac   dE  = geom->dimEmbed;
7449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
7459566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS));
7469566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(field, cellIS, NULL, &maxDegree));
7479566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointIS, &points));
7480145028aSToby Isaac   numFaces = geom->numCells;
7490145028aSToby Isaac   Nq       = geom->numPoints;
750f15274beSMatthew Knepley   /* First, set local faces and flip normals so that they are outward for the first supporting cell */
751f15274beSMatthew Knepley   for (p = 0; p < numFaces; p++) {
752f15274beSMatthew Knepley     PetscInt        point = points[p];
753f15274beSMatthew Knepley     PetscInt        suppSize, s, coneSize, c, numChildren;
754f15274beSMatthew Knepley     const PetscInt *supp, *cone, *ornt;
755f15274beSMatthew Knepley 
7569566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
75728b400f6SJacob Faibussowitsch     PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
7589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, point, &suppSize));
75963a3b9bcSJacob Faibussowitsch     PetscCheck(suppSize <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, suppSize);
760f15274beSMatthew Knepley     if (!suppSize) continue;
7619566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, point, &supp));
762f15274beSMatthew Knepley     for (s = 0; s < suppSize; ++s) {
7639566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, supp[s], &coneSize));
7649566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, supp[s], &cone));
7659371c9d4SSatish Balay       for (c = 0; c < coneSize; ++c)
7669371c9d4SSatish Balay         if (cone[c] == point) break;
76763a3b9bcSJacob 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]);
768f15274beSMatthew Knepley       geom->face[p][s] = c;
769f15274beSMatthew Knepley     }
7709566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeOrientation(dm, supp[0], &ornt));
771f15274beSMatthew Knepley     if (ornt[geom->face[p][0]] < 0) {
772f15274beSMatthew Knepley       PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d;
773f15274beSMatthew Knepley 
7749371c9d4SSatish Balay       for (q = 0; q < Np; ++q)
7759371c9d4SSatish Balay         for (d = 0; d < dE; ++d) geom->n[(p * Np + q) * dE + d] = -geom->n[(p * Np + q) * dE + d];
776f15274beSMatthew Knepley     }
777f15274beSMatthew Knepley   }
778b7260050SToby Isaac   if (maxDegree <= 1) {
7790145028aSToby Isaac     PetscInt     numCells, offset, *cells;
7800145028aSToby Isaac     PetscFEGeom *cellGeom;
7810145028aSToby Isaac     IS           suppIS;
7820145028aSToby Isaac 
7830145028aSToby Isaac     for (p = 0, numCells = 0; p < numFaces; p++) {
7840145028aSToby Isaac       PetscInt point = points[p];
7850145028aSToby Isaac       PetscInt numSupp, numChildren;
7860145028aSToby Isaac 
7879566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
78828b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
7899566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
79063a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
7910145028aSToby Isaac       numCells += numSupp;
7920145028aSToby Isaac     }
7939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numCells, &cells));
7940145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
7950145028aSToby Isaac       PetscInt        point = points[p];
7960145028aSToby Isaac       PetscInt        numSupp, s;
7970145028aSToby Isaac       const PetscInt *supp;
7980145028aSToby Isaac 
7999566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8009566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
801ad540459SPierre Jolivet       for (s = 0; s < numSupp; s++, offset++) cells[offset] = supp[s];
8020145028aSToby Isaac     }
8039566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
8049566063dSJacob Faibussowitsch     PetscCall(DMFieldCreateFEGeom(field, suppIS, quad, PETSC_FALSE, &cellGeom));
8050145028aSToby Isaac     for (p = 0, offset = 0; p < numFaces; p++) {
8060145028aSToby Isaac       PetscInt        point = points[p];
8070145028aSToby Isaac       PetscInt        numSupp, s, q;
8080145028aSToby Isaac       const PetscInt *supp;
8090145028aSToby Isaac 
8109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
8119566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
8120145028aSToby Isaac       for (s = 0; s < numSupp; s++, offset++) {
8130145028aSToby Isaac         for (q = 0; q < Nq * dE * dE; q++) {
81441418987SMatthew G. Knepley           geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
8150145028aSToby Isaac           geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
8160145028aSToby Isaac         }
81741418987SMatthew G. Knepley         for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
8180145028aSToby Isaac       }
8190145028aSToby Isaac     }
8209566063dSJacob Faibussowitsch     PetscCall(PetscFEGeomDestroy(&cellGeom));
8219566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&suppIS));
8229566063dSJacob Faibussowitsch     PetscCall(PetscFree(cells));
8230145028aSToby Isaac   } else {
8246858538eSMatthew G. Knepley     DMField_DS    *dsfield = (DMField_DS *)field->data;
8250145028aSToby Isaac     PetscObject    faceDisc, cellDisc;
8260145028aSToby Isaac     PetscClassId   faceId, cellId;
8270145028aSToby Isaac     PetscDualSpace dsp;
8280145028aSToby Isaac     DM             K;
829ba2698f1SMatthew G. Knepley     DMPolytopeType ct;
8300145028aSToby Isaac     PetscInt(*co)[2][3];
8310145028aSToby Isaac     PetscInt        coneSize;
8320145028aSToby Isaac     PetscInt      **counts;
8330145028aSToby Isaac     PetscInt        f, i, o, q, s;
8340145028aSToby Isaac     const PetscInt *coneK;
835ba2698f1SMatthew G. Knepley     PetscInt        eStart, minOrient, maxOrient, numOrient;
8360145028aSToby Isaac     PetscInt       *orients;
8370145028aSToby Isaac     PetscReal     **orientPoints;
8380145028aSToby Isaac     PetscReal      *cellPoints;
8390145028aSToby Isaac     PetscReal      *dummyWeights;
8400145028aSToby Isaac     PetscQuadrature cellQuad = NULL;
8410145028aSToby Isaac 
8426858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 1, dsfield->disc, &faceDisc));
8436858538eSMatthew G. Knepley     PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc));
8449566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(faceDisc, &faceId));
8459566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(cellDisc, &cellId));
8461dca8a05SBarry Smith     PetscCheck(faceId == PETSCFE_CLASSID && cellId == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported");
8479566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace((PetscFE)cellDisc, &dsp));
8489566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(dsp, &K));
8499566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(K, 1, &eStart, NULL));
8509566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(K, eStart, &ct));
8519566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(K, 0, &coneSize));
8529566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
8539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numFaces, &co, coneSize, &counts));
8549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dE * Nq, &cellPoints));
8559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nq, &dummyWeights));
8569566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad));
8579566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights));
8580145028aSToby Isaac     minOrient = PETSC_MAX_INT;
8590145028aSToby Isaac     maxOrient = PETSC_MIN_INT;
8600145028aSToby Isaac     for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */
8610145028aSToby Isaac       PetscInt        point = points[p];
8620145028aSToby Isaac       PetscInt        numSupp, numChildren;
8630145028aSToby Isaac       const PetscInt *supp;
8640145028aSToby Isaac 
8659566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL));
86628b400f6SJacob Faibussowitsch       PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
8679566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm, point, &numSupp));
86863a3b9bcSJacob Faibussowitsch       PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp);
8699566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupport(dm, point, &supp));
8700145028aSToby Isaac       for (s = 0; s < numSupp; s++) {
8710145028aSToby Isaac         PetscInt        cell = supp[s];
8720145028aSToby Isaac         PetscInt        numCone;
8730145028aSToby Isaac         const PetscInt *cone, *orient;
8740145028aSToby Isaac 
8759566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(dm, cell, &numCone));
87608401ef6SPierre Jolivet         PetscCheck(numCone == coneSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element");
8779566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(dm, cell, &cone));
8789566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeOrientation(dm, cell, &orient));
8790145028aSToby Isaac         for (f = 0; f < coneSize; f++) {
8800145028aSToby Isaac           if (cone[f] == point) break;
8810145028aSToby Isaac         }
8820145028aSToby Isaac         co[p][s][0] = f;
8830145028aSToby Isaac         co[p][s][1] = orient[f];
8840145028aSToby Isaac         co[p][s][2] = cell;
8850145028aSToby Isaac         minOrient   = PetscMin(minOrient, orient[f]);
886cd93b0e1SMatthew G. Knepley         maxOrient   = PetscMax(maxOrient, orient[f]);
8870145028aSToby Isaac       }
8880145028aSToby Isaac       for (; s < 2; s++) {
8890145028aSToby Isaac         co[p][s][0] = -1;
8900145028aSToby Isaac         co[p][s][1] = -1;
8910145028aSToby Isaac         co[p][s][2] = -1;
8920145028aSToby Isaac       }
8930145028aSToby Isaac     }
8940145028aSToby Isaac     numOrient = maxOrient + 1 - minOrient;
8959566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(K, 0, &coneK));
8960145028aSToby Isaac     /* count all (face,orientation) doubles that appear */
8979566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(numOrient, &orients, numOrient, &orientPoints));
8989566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscCalloc1(numOrient + 1, &counts[f]));
8990145028aSToby Isaac     for (p = 0; p < numFaces; p++) {
9000145028aSToby Isaac       for (s = 0; s < 2; s++) {
9010145028aSToby Isaac         if (co[p][s][0] >= 0) {
9020145028aSToby Isaac           counts[co[p][s][0]][co[p][s][1] - minOrient]++;
9030145028aSToby Isaac           orients[co[p][s][1] - minOrient]++;
9040145028aSToby Isaac         }
9050145028aSToby Isaac       }
9060145028aSToby Isaac     }
9070145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
9080145028aSToby Isaac       if (orients[o]) {
9090145028aSToby Isaac         PetscInt orient = o + minOrient;
9100145028aSToby Isaac         PetscInt q;
9110145028aSToby Isaac 
9129566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Nq * dim, &orientPoints[o]));
9130145028aSToby Isaac         /* rotate the quadrature points appropriately */
914ba2698f1SMatthew G. Knepley         switch (ct) {
915ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_POINT: break;
916ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_SEGMENT:
9170145028aSToby Isaac           if (orient == -2 || orient == 1) {
918ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = -geom->xi[q];
9190145028aSToby Isaac           } else {
920ad540459SPierre Jolivet             for (q = 0; q < Nq; q++) orientPoints[o][q] = geom->xi[q];
9210145028aSToby Isaac           }
9220145028aSToby Isaac           break;
923ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_TRIANGLE:
9240145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9250145028aSToby Isaac             PetscReal lambda[3];
9260145028aSToby Isaac             PetscReal lambdao[3];
9270145028aSToby Isaac 
9280145028aSToby Isaac             /* convert to barycentric */
9290145028aSToby Isaac             lambda[0] = -(geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.;
9300145028aSToby Isaac             lambda[1] = (geom->xi[2 * q] + 1.) / 2.;
9310145028aSToby Isaac             lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.;
9320145028aSToby Isaac             if (orient >= 0) {
933ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(orient + i) % 3];
9340145028aSToby Isaac             } else {
935ad540459SPierre Jolivet               for (i = 0; i < 3; i++) lambdao[i] = lambda[(-(orient + i) + 3) % 3];
9360145028aSToby Isaac             }
9370145028aSToby Isaac             /* convert to coordinates */
9380145028aSToby Isaac             orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1];
9390145028aSToby Isaac             orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2];
9400145028aSToby Isaac           }
9410145028aSToby Isaac           break;
942ba2698f1SMatthew G. Knepley         case DM_POLYTOPE_QUADRILATERAL:
9430145028aSToby Isaac           for (q = 0; q < Nq; q++) {
9440145028aSToby Isaac             PetscReal xi[2], xio[2];
9450145028aSToby Isaac             PetscInt  oabs = (orient >= 0) ? orient : -(orient + 1);
9460145028aSToby Isaac 
9470145028aSToby Isaac             xi[0] = geom->xi[2 * q];
9480145028aSToby Isaac             xi[1] = geom->xi[2 * q + 1];
9490145028aSToby Isaac             switch (oabs) {
9500145028aSToby Isaac             case 1:
9510145028aSToby Isaac               xio[0] = xi[1];
9520145028aSToby Isaac               xio[1] = -xi[0];
9530145028aSToby Isaac               break;
9549371c9d4SSatish Balay             case 2: xio[0] = -xi[0]; xio[1] = -xi[1];
9559371c9d4SSatish Balay             case 3: xio[0] = -xi[1]; xio[1] = xi[0];
956aa1371cbSToby Isaac             case 0:
957aa1371cbSToby Isaac             default:
958aa1371cbSToby Isaac               xio[0] = xi[0];
959aa1371cbSToby Isaac               xio[1] = xi[1];
960aa1371cbSToby Isaac               break;
9610145028aSToby Isaac             }
962ad540459SPierre Jolivet             if (orient < 0) xio[0] = -xio[0];
9630145028aSToby Isaac             orientPoints[o][2 * q + 0] = xio[0];
9640145028aSToby Isaac             orientPoints[o][2 * q + 1] = xio[1];
9650145028aSToby Isaac           }
9660145028aSToby Isaac           break;
96798921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s not yet supported", DMPolytopeTypes[ct]);
9680145028aSToby Isaac         }
9690145028aSToby Isaac       }
9700145028aSToby Isaac     }
9710145028aSToby Isaac     for (f = 0; f < coneSize; f++) {
9720145028aSToby Isaac       PetscInt  face = coneK[f];
9734d1a973fSToby Isaac       PetscReal v0[3];
9742d89661fSMatthew G. Knepley       PetscReal J[9], detJ;
9750145028aSToby Isaac       PetscInt  numCells, offset;
9760145028aSToby Isaac       PetscInt *cells;
9770145028aSToby Isaac       IS        suppIS;
9780145028aSToby Isaac 
9799566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ));
9800145028aSToby Isaac       for (o = 0; o <= numOrient; o++) {
9810145028aSToby Isaac         PetscFEGeom *cellGeom;
9820145028aSToby Isaac 
9830145028aSToby Isaac         if (!counts[f][o]) continue;
9840145028aSToby Isaac         /* If this (face,orientation) double appears,
9850145028aSToby Isaac          * convert the face quadrature points into volume quadrature points */
9860145028aSToby Isaac         for (q = 0; q < Nq; q++) {
9870145028aSToby Isaac           PetscReal xi0[3] = {-1., -1., -1.};
9880145028aSToby Isaac 
9890145028aSToby Isaac           CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]);
9900145028aSToby Isaac         }
9910145028aSToby Isaac         for (p = 0, numCells = 0; p < numFaces; p++) {
9920145028aSToby Isaac           for (s = 0; s < 2; s++) {
9930145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++;
9940145028aSToby Isaac           }
9950145028aSToby Isaac         }
9969566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(numCells, &cells));
9970145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
9980145028aSToby Isaac           for (s = 0; s < 2; s++) {
999ad540459SPierre Jolivet             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) cells[offset++] = co[p][s][2];
10000145028aSToby Isaac           }
10010145028aSToby Isaac         }
10029566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS));
10039566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FALSE, &cellGeom));
10040145028aSToby Isaac         for (p = 0, offset = 0; p < numFaces; p++) {
10050145028aSToby Isaac           for (s = 0; s < 2; s++) {
10060145028aSToby Isaac             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
10070145028aSToby Isaac               for (q = 0; q < Nq * dE * dE; q++) {
100841418987SMatthew G. Knepley                 geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
10090145028aSToby Isaac                 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
10100145028aSToby Isaac               }
101141418987SMatthew G. Knepley               for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
10120145028aSToby Isaac               offset++;
10130145028aSToby Isaac             }
10140145028aSToby Isaac           }
10150145028aSToby Isaac         }
10169566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&cellGeom));
10179566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&suppIS));
10189566063dSJacob Faibussowitsch         PetscCall(PetscFree(cells));
10190145028aSToby Isaac       }
10200145028aSToby Isaac     }
10210145028aSToby Isaac     for (o = 0; o < numOrient; o++) {
102248a46eb9SPierre Jolivet       if (orients[o]) PetscCall(PetscFree(orientPoints[o]));
10230145028aSToby Isaac     }
10249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(orients, orientPoints));
10259566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&cellQuad));
10269566063dSJacob Faibussowitsch     for (f = 0; f < coneSize; f++) PetscCall(PetscFree(counts[f]));
10279566063dSJacob Faibussowitsch     PetscCall(PetscFree2(co, counts));
10280145028aSToby Isaac   }
10299566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointIS, &points));
10309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
10310145028aSToby Isaac   PetscFunctionReturn(0);
10320145028aSToby Isaac }
10330145028aSToby Isaac 
10349371c9d4SSatish Balay static PetscErrorCode DMFieldInitialize_DS(DMField field) {
103551a454edSToby Isaac   PetscFunctionBegin;
103651a454edSToby Isaac   field->ops->destroy                 = DMFieldDestroy_DS;
103751a454edSToby Isaac   field->ops->evaluate                = DMFieldEvaluate_DS;
103851a454edSToby Isaac   field->ops->evaluateFE              = DMFieldEvaluateFE_DS;
1039805e7170SToby Isaac   field->ops->evaluateFV              = DMFieldEvaluateFV_DS;
1040b7260050SToby Isaac   field->ops->getDegree               = DMFieldGetDegree_DS;
104151a454edSToby Isaac   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS;
104251a454edSToby Isaac   field->ops->view                    = DMFieldView_DS;
10430145028aSToby Isaac   field->ops->computeFaceData         = DMFieldComputeFaceData_DS;
104451a454edSToby Isaac   PetscFunctionReturn(0);
104551a454edSToby Isaac }
104651a454edSToby Isaac 
10479371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field) {
104851a454edSToby Isaac   DMField_DS *dsfield;
104951a454edSToby Isaac 
105051a454edSToby Isaac   PetscFunctionBegin;
1051*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dsfield));
105251a454edSToby Isaac   field->data = dsfield;
10539566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_DS(field));
105451a454edSToby Isaac   PetscFunctionReturn(0);
105551a454edSToby Isaac }
105651a454edSToby Isaac 
10579371c9d4SSatish Balay PetscErrorCode DMFieldCreateDSWithDG(DM dm, DM dmDG, PetscInt fieldNum, Vec vec, Vec vecDG, DMField *field) {
105851a454edSToby Isaac   DMField      b;
105951a454edSToby Isaac   DMField_DS  *dsfield;
10606858538eSMatthew G. Knepley   PetscObject  disc = NULL, discDG = NULL;
10616858538eSMatthew G. Knepley   PetscSection section;
106251a454edSToby Isaac   PetscBool    isContainer   = PETSC_FALSE;
106351a454edSToby Isaac   PetscClassId id            = -1;
106451a454edSToby Isaac   PetscInt     numComponents = -1, dsNumFields;
106551a454edSToby Isaac 
106651a454edSToby Isaac   PetscFunctionBegin;
10679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
10689566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetFieldComponents(section, fieldNum, &numComponents));
10699566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &dsNumFields));
10709566063dSJacob Faibussowitsch   if (dsNumFields) PetscCall(DMGetField(dm, fieldNum, NULL, &disc));
10716858538eSMatthew G. Knepley   if (dsNumFields && dmDG) {
10726858538eSMatthew G. Knepley     PetscCall(DMGetField(dmDG, fieldNum, NULL, &discDG));
10736858538eSMatthew G. Knepley     PetscCall(PetscObjectReference(discDG));
10746858538eSMatthew G. Knepley   }
107551a454edSToby Isaac   if (disc) {
10769566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(disc, &id));
107751a454edSToby Isaac     isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE;
107851a454edSToby Isaac   }
107951a454edSToby Isaac   if (!disc || isContainer) {
108051a454edSToby Isaac     MPI_Comm       comm = PetscObjectComm((PetscObject)dm);
108151a454edSToby Isaac     PetscFE        fe;
10822df84da0SMatthew G. Knepley     DMPolytopeType ct, locct = DM_POLYTOPE_UNKNOWN;
10832df84da0SMatthew G. Knepley     PetscInt       dim, cStart, cEnd, cellHeight;
108451a454edSToby Isaac 
10859566063dSJacob Faibussowitsch     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
10869566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
10879566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
10889566063dSJacob Faibussowitsch     if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &locct));
10899566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&locct, &ct, 1, MPI_INT, MPI_MIN, comm));
10909566063dSJacob Faibussowitsch     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, numComponents, ct, 1, PETSC_DETERMINE, &fe));
10919566063dSJacob Faibussowitsch     PetscCall(PetscFEViewFromOptions(fe, NULL, "-field_fe_view"));
109251a454edSToby Isaac     disc = (PetscObject)fe;
10931baa6e33SBarry Smith   } else PetscCall(PetscObjectReference(disc));
10949566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(disc, &id));
10956858538eSMatthew G. Knepley   if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)disc, &numComponents));
10966858538eSMatthew G. Knepley   else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine number of discretization components");
10979566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm, numComponents, DMFIELD_VERTEX, &b));
10989566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b, DMFIELDDS));
109951a454edSToby Isaac   dsfield           = (DMField_DS *)b->data;
110051a454edSToby Isaac   dsfield->fieldNum = fieldNum;
11019566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dsfield->height));
1102f99c8401SToby Isaac   dsfield->height++;
11039566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(dsfield->height, &dsfield->disc));
110451a454edSToby Isaac   dsfield->disc[0] = disc;
11059566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)vec));
110651a454edSToby Isaac   dsfield->vec = vec;
11076858538eSMatthew G. Knepley   if (dmDG) {
11086858538eSMatthew G. Knepley     dsfield->dmDG = dmDG;
11096858538eSMatthew G. Knepley     PetscCall(PetscCalloc1(dsfield->height, &dsfield->discDG));
11106858538eSMatthew G. Knepley     dsfield->discDG[0] = discDG;
11116858538eSMatthew G. Knepley     PetscCall(PetscObjectReference((PetscObject)vecDG));
11126858538eSMatthew G. Knepley     dsfield->vecDG = vecDG;
11136858538eSMatthew G. Knepley   }
111451a454edSToby Isaac   *field = b;
11156858538eSMatthew G. Knepley   PetscFunctionReturn(0);
11166858538eSMatthew G. Knepley }
111751a454edSToby Isaac 
11189371c9d4SSatish Balay PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec, DMField *field) {
11196858538eSMatthew G. Knepley   PetscFunctionBegin;
11206858538eSMatthew G. Knepley   PetscCall(DMFieldCreateDSWithDG(dm, NULL, fieldNum, vec, NULL, field));
112151a454edSToby Isaac   PetscFunctionReturn(0);
112251a454edSToby Isaac }
1123