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; 39*9f196a02SMartin Diehl PetscBool isascii; 4051a454edSToby Isaac 4151a454edSToby Isaac PetscFunctionBegin; 42*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 4351a454edSToby Isaac disc = dsfield->disc[0]; 44*9f196a02SMartin Diehl if (isascii) { 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)); 139835f2295SStefano Zampini PetscCall(DMPlexVecRestoreClosure(fdm, s, dsfield->vec, cell, Nc, 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, §ion)); 1659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(section, dsfield->fieldNum, §ion)); 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; 250835f2295SStefano Zampini PetscMPIInt idim; 251a6216909SToby Isaac 252a6216909SToby Isaac PetscFunctionBegin; 253a6216909SToby Isaac nc = field->numComponents; 2549566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(field->dm, §ion)); 2556858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(cellDisc, &discID)); 25708401ef6SPierre Jolivet PetscCheck(discID == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization type not supported"); 258a6216909SToby Isaac cellFE = (PetscFE)cellDisc; 2599566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(cellFE, &feDim)); 2609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(field->dm, &dim)); 2619566063dSJacob Faibussowitsch PetscCall(DMGetDimension(field->dm, &dimR)); 2629566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells)); 264ad540459SPierre 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); 2659566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(cellSF, &cellDegrees)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(cellSF, &cellDegrees)); 267a6216909SToby Isaac for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) { 268a6216909SToby Isaac gatherMax = PetscMax(gatherMax, cellDegrees[c]); 269a6216909SToby Isaac gatherSize += cellDegrees[c]; 270a6216909SToby Isaac } 2719566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(gatherSize * dim, &cellPoints, gatherMax * dim, &coordsReal, gatherMax * dimR, &coordsRef)); 2729566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(gatherMax * dimR, &v, gatherMax * dimR * dimR, &J, gatherMax * dimR * dimR, &invJ, gatherMax, &detJ)); 27332603206SJames Wright if (datatype == PETSC_SCALAR) PetscCall(PetscMalloc3((B ? (size_t)nc * gatherSize : 0), &cellBs, (D ? (size_t)nc * dim * gatherSize : 0), &cellDs, (H ? (size_t)nc * dim * dim * gatherSize : 0), &cellHs)); 27432603206SJames Wright else PetscCall(PetscMalloc3((B ? (size_t)nc * gatherSize : 0), &cellBr, (D ? (size_t)nc * dim * gatherSize : 0), &cellDr, (H ? (size_t)nc * dim * dim * gatherSize : 0), &cellHr)); 275a6216909SToby Isaac 276835f2295SStefano Zampini PetscCall(PetscMPIIntCast(dim, &idim)); 277835f2295SStefano Zampini PetscCallMPI(MPI_Type_contiguous(idim, MPIU_SCALAR, &pointType)); 2789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&pointType)); 2799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(points, &pointsArray)); 2809566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints)); 2819566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints)); 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(points, &pointsArray)); 283a6216909SToby Isaac for (c = 0, offset = 0; c < numCells; c++) { 284a6216909SToby Isaac PetscInt nq = cellDegrees[c], p; 285a6216909SToby Isaac 286a6216909SToby Isaac if (nq) { 287ef0bb6c7SMatthew G. Knepley PetscInt K = H ? 2 : (D ? 1 : (B ? 0 : -1)); 288ef0bb6c7SMatthew G. Knepley PetscTabulation T; 2896858538eSMatthew G. Knepley PetscQuadrature quad; 2906858538eSMatthew G. Knepley const PetscScalar *array; 291a6216909SToby Isaac PetscScalar *elem = NULL; 2924cbe5319SToby Isaac PetscReal *quadPoints; 2936858538eSMatthew G. Knepley PetscBool isDG; 2946858538eSMatthew G. Knepley PetscInt closureSize, d, e, f, g; 295a6216909SToby Isaac 2964d1a973fSToby Isaac for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]); 2979566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef)); 2989566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(cellFE, 1, nq, coordsRef, K, &T)); 2999566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad)); 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dimR * nq, &quadPoints)); 3014cbe5319SToby Isaac for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p]; 3029566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL)); 3039566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ)); 3049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 3056858538eSMatthew G. Knepley PetscCall(DMFieldGetClosure_Internal(field, c, &isDG, &closureSize, &array, &elem)); 306a6216909SToby Isaac if (B) { 307a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 308a6216909SToby Isaac PetscScalar *cB = &cellBs[nc * offset]; 309a6216909SToby Isaac 310ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, (PetscScalar)); 311a6216909SToby Isaac } else { 312a6216909SToby Isaac PetscReal *cB = &cellBr[nc * offset]; 313a6216909SToby Isaac 314ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cB, T->T[0], elem, nq, feDim, nc, PetscRealPart); 315a6216909SToby Isaac } 316a6216909SToby Isaac } 317a6216909SToby Isaac if (D) { 318a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 319a6216909SToby Isaac PetscScalar *cD = &cellDs[nc * dim * offset]; 320a6216909SToby Isaac 321ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), (PetscScalar)); 3224cbe5319SToby Isaac for (p = 0; p < nq; p++) { 3234cbe5319SToby Isaac for (g = 0; g < nc; g++) { 3244d1a973fSToby Isaac PetscScalar vs[3]; 3254d1a973fSToby Isaac 3264cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 3274d1a973fSToby Isaac vs[d] = 0.; 328ad540459SPierre Jolivet for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 3294cbe5319SToby Isaac } 330ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = vs[d]; 3314cbe5319SToby Isaac } 3324cbe5319SToby Isaac } 333a6216909SToby Isaac } else { 334a6216909SToby Isaac PetscReal *cD = &cellDr[nc * dim * offset]; 335a6216909SToby Isaac 336ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cD, T->T[1], elem, nq, feDim, (nc * dim), PetscRealPart); 3374cbe5319SToby Isaac for (p = 0; p < nq; p++) { 3384cbe5319SToby Isaac for (g = 0; g < nc; g++) { 3394cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 3404cbe5319SToby Isaac v[d] = 0.; 341ad540459SPierre Jolivet for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e]; 3424cbe5319SToby Isaac } 343ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cD[(nc * p + g) * dimR + d] = v[d]; 3444cbe5319SToby Isaac } 3454cbe5319SToby Isaac } 346a6216909SToby Isaac } 347a6216909SToby Isaac } 348a6216909SToby Isaac if (H) { 349a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 350a6216909SToby Isaac PetscScalar *cH = &cellHs[nc * dim * dim * offset]; 351a6216909SToby Isaac 352ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), (PetscScalar)); 3534cbe5319SToby Isaac for (p = 0; p < nq; p++) { 3544cbe5319SToby Isaac for (g = 0; g < nc * dimR; g++) { 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 * dimR * p + g) * dimR + e]; 3604cbe5319SToby Isaac } 361ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = vs[d]; 3624cbe5319SToby Isaac } 3634cbe5319SToby Isaac for (g = 0; g < nc; g++) { 3644cbe5319SToby Isaac for (f = 0; f < dimR; f++) { 3654d1a973fSToby Isaac PetscScalar vs[3]; 3664d1a973fSToby Isaac 3674cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 3684d1a973fSToby Isaac vs[d] = 0.; 369ad540459SPierre Jolivet for (e = 0; e < dimR; e++) vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 3704cbe5319SToby Isaac } 371ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d]; 3724cbe5319SToby Isaac } 3734cbe5319SToby Isaac } 3744cbe5319SToby Isaac } 375a6216909SToby Isaac } else { 376a6216909SToby Isaac PetscReal *cH = &cellHr[nc * dim * dim * offset]; 377a6216909SToby Isaac 378ef0bb6c7SMatthew G. Knepley DMFieldDSdot(cH, T->T[2], elem, nq, feDim, (nc * dim * dim), PetscRealPart); 3794cbe5319SToby Isaac for (p = 0; p < nq; p++) { 3804cbe5319SToby Isaac for (g = 0; g < nc * dimR; g++) { 3814cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 3824cbe5319SToby Isaac v[d] = 0.; 383ad540459SPierre Jolivet for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e]; 3844cbe5319SToby Isaac } 385ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cH[(nc * dimR * p + g) * dimR + d] = v[d]; 3864cbe5319SToby Isaac } 3874cbe5319SToby Isaac for (g = 0; g < nc; g++) { 3884cbe5319SToby Isaac for (f = 0; f < dimR; f++) { 3894cbe5319SToby Isaac for (d = 0; d < dimR; d++) { 3904cbe5319SToby Isaac v[d] = 0.; 391ad540459SPierre Jolivet for (e = 0; e < dimR; e++) v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f]; 3924cbe5319SToby Isaac } 393ad540459SPierre Jolivet for (d = 0; d < dimR; d++) cH[((nc * p + g) * dimR + d) * dimR + f] = v[d]; 3944cbe5319SToby Isaac } 3954cbe5319SToby Isaac } 3964cbe5319SToby Isaac } 397a6216909SToby Isaac } 398a6216909SToby Isaac } 3996858538eSMatthew G. Knepley PetscCall(DMFieldRestoreClosure_Internal(field, c, &isDG, &closureSize, &array, &elem)); 4009566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&T)); 401a6216909SToby Isaac } 402a6216909SToby Isaac offset += nq; 403a6216909SToby Isaac } 404a6216909SToby Isaac { 405a6216909SToby Isaac MPI_Datatype origtype; 4064cbe5319SToby Isaac 407a6216909SToby Isaac if (datatype == PETSC_SCALAR) { 408a6216909SToby Isaac origtype = MPIU_SCALAR; 409a6216909SToby Isaac } else { 410a6216909SToby Isaac origtype = MPIU_REAL; 411a6216909SToby Isaac } 412a6216909SToby Isaac if (B) { 413a6216909SToby Isaac MPI_Datatype Btype; 414a6216909SToby Isaac 415835f2295SStefano Zampini PetscCall(PetscMPIIntCast(nc, &idim)); 416835f2295SStefano Zampini PetscCallMPI(MPI_Type_contiguous(idim, origtype, &Btype)); 4179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&Btype)); 4189566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B)); 4199566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(cellSF, Btype, (datatype == PETSC_SCALAR) ? (void *)cellBs : (void *)cellBr, B)); 4209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&Btype)); 421a6216909SToby Isaac } 422a6216909SToby Isaac if (D) { 423a6216909SToby Isaac MPI_Datatype Dtype; 424a6216909SToby Isaac 425835f2295SStefano Zampini PetscCall(PetscMPIIntCast(nc * dim, &idim)); 426835f2295SStefano Zampini PetscCallMPI(MPI_Type_contiguous(idim, origtype, &Dtype)); 4279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&Dtype)); 4289566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D)); 4299566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(cellSF, Dtype, (datatype == PETSC_SCALAR) ? (void *)cellDs : (void *)cellDr, D)); 4309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&Dtype)); 431a6216909SToby Isaac } 432a6216909SToby Isaac if (H) { 433a6216909SToby Isaac MPI_Datatype Htype; 434a6216909SToby Isaac 435835f2295SStefano Zampini PetscCall(PetscMPIIntCast(nc * dim * dim, &idim)); 436835f2295SStefano Zampini PetscCallMPI(MPI_Type_contiguous(idim, origtype, &Htype)); 4379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&Htype)); 4389566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H)); 4399566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(cellSF, Htype, (datatype == PETSC_SCALAR) ? (void *)cellHs : (void *)cellHr, H)); 4409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&Htype)); 441a6216909SToby Isaac } 442a6216909SToby Isaac } 4439566063dSJacob Faibussowitsch PetscCall(PetscFree4(v, J, invJ, detJ)); 4449566063dSJacob Faibussowitsch PetscCall(PetscFree3(cellBr, cellDr, cellHr)); 4459566063dSJacob Faibussowitsch PetscCall(PetscFree3(cellBs, cellDs, cellHs)); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree3(cellPoints, coordsReal, coordsRef)); 4479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&pointType)); 4489566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&cellSF)); 4493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 450a6216909SToby Isaac } 451a6216909SToby Isaac 452d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) 453d71ae5a4SJacob Faibussowitsch { 454805e7170SToby Isaac DMField_DS *dsfield = (DMField_DS *)field->data; 455805e7170SToby Isaac PetscInt h, imin; 456805e7170SToby Isaac PetscInt dim; 457805e7170SToby Isaac PetscClassId id; 458805e7170SToby Isaac PetscQuadrature quad = NULL; 459b7260050SToby Isaac PetscInt maxDegree; 460805e7170SToby Isaac PetscFEGeom *geom; 461805e7170SToby Isaac PetscInt Nq, Nc, dimC, qNc, N; 462805e7170SToby Isaac PetscInt numPoints; 463805e7170SToby Isaac void *qB = NULL, *qD = NULL, *qH = NULL; 464805e7170SToby Isaac const PetscReal *weights; 465805e7170SToby Isaac MPI_Datatype mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL; 466805e7170SToby Isaac PetscObject disc; 467805e7170SToby Isaac DMField coordField; 468805e7170SToby Isaac 469805e7170SToby Isaac PetscFunctionBegin; 470805e7170SToby Isaac Nc = field->numComponents; 4719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(field->dm, &dimC)); 4729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(field->dm, &dim)); 4739566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numPoints)); 4749566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(pointIS, &imin, NULL)); 475805e7170SToby Isaac for (h = 0; h < dsfield->height; h++) { 476805e7170SToby Isaac PetscInt hEnd; 477805e7170SToby Isaac 4789566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd)); 479805e7170SToby Isaac if (imin < hEnd) break; 480805e7170SToby Isaac } 481805e7170SToby Isaac dim -= h; 4826858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc)); 4839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 48408401ef6SPierre Jolivet PetscCheck(id == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported"); 4859566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateField(field->dm, &coordField)); 4869566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree)); 48748a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad)); 4889566063dSJacob Faibussowitsch if (!quad) PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, &quad)); 48928b400f6SJacob Faibussowitsch PetscCheck(quad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages"); 490ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FEGEOM_BASIC, &geom)); 4919566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights)); 49208401ef6SPierre Jolivet PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components"); 493805e7170SToby Isaac N = numPoints * Nq * Nc; 4949566063dSJacob Faibussowitsch if (B) PetscCall(DMGetWorkArray(field->dm, N, mpitype, &qB)); 4959566063dSJacob Faibussowitsch if (D) PetscCall(DMGetWorkArray(field->dm, N * dimC, mpitype, &qD)); 4969566063dSJacob Faibussowitsch if (H) PetscCall(DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH)); 4979566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field, pointIS, quad, type, qB, qD, qH)); 498805e7170SToby Isaac if (B) { 499805e7170SToby Isaac PetscInt i, j, k; 500805e7170SToby Isaac 501805e7170SToby Isaac if (type == PETSC_SCALAR) { 502805e7170SToby Isaac PetscScalar *sB = (PetscScalar *)B; 503805e7170SToby Isaac PetscScalar *sqB = (PetscScalar *)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++) sB[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++) sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[(i * Nq + k) * Nc + j]; 512805e7170SToby Isaac } 51362bd480fSMatthew G. Knepley for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol; 514805e7170SToby Isaac } 515805e7170SToby Isaac } else { 516805e7170SToby Isaac PetscReal *rB = (PetscReal *)B; 517805e7170SToby Isaac PetscReal *rqB = (PetscReal *)qB; 518805e7170SToby Isaac 519805e7170SToby Isaac for (i = 0; i < numPoints; i++) { 520805e7170SToby Isaac PetscReal vol = 0.; 521805e7170SToby Isaac 522ad540459SPierre Jolivet for (j = 0; j < Nc; j++) rB[i * Nc + j] = 0.; 523805e7170SToby Isaac for (k = 0; k < Nq; k++) { 524805e7170SToby Isaac vol += geom->detJ[i * Nq + k] * weights[k]; 525ad540459SPierre Jolivet for (j = 0; j < Nc; j++) rB[i * Nc + j] += weights[k] * rqB[(i * Nq + k) * Nc + j]; 526805e7170SToby Isaac } 527805e7170SToby Isaac for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol; 528805e7170SToby Isaac } 529805e7170SToby Isaac } 530805e7170SToby Isaac } 531805e7170SToby Isaac if (D) { 532805e7170SToby Isaac PetscInt i, j, k, l, m; 533805e7170SToby Isaac 534805e7170SToby Isaac if (type == PETSC_SCALAR) { 535805e7170SToby Isaac PetscScalar *sD = (PetscScalar *)D; 536805e7170SToby Isaac PetscScalar *sqD = (PetscScalar *)qD; 537805e7170SToby Isaac 538805e7170SToby Isaac for (i = 0; i < numPoints; i++) { 539805e7170SToby Isaac PetscReal vol = 0.; 540805e7170SToby Isaac 541ad540459SPierre Jolivet for (j = 0; j < Nc * dimC; j++) sD[i * Nc * dimC + j] = 0.; 542805e7170SToby Isaac for (k = 0; k < Nq; k++) { 543805e7170SToby Isaac vol += geom->detJ[i * Nq + k] * weights[k]; 544805e7170SToby Isaac for (j = 0; j < Nc; j++) { 545805e7170SToby Isaac PetscScalar pD[3] = {0., 0., 0.}; 546805e7170SToby Isaac 547805e7170SToby Isaac for (l = 0; l < dimC; l++) { 548ad540459SPierre 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]; 549805e7170SToby Isaac } 550ad540459SPierre Jolivet for (l = 0; l < dimC; l++) sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l]; 551805e7170SToby Isaac } 552805e7170SToby Isaac } 553805e7170SToby Isaac for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol; 554805e7170SToby Isaac } 555805e7170SToby Isaac } else { 556805e7170SToby Isaac PetscReal *rD = (PetscReal *)D; 557805e7170SToby Isaac PetscReal *rqD = (PetscReal *)qD; 558805e7170SToby Isaac 559805e7170SToby Isaac for (i = 0; i < numPoints; i++) { 560805e7170SToby Isaac PetscReal vol = 0.; 561805e7170SToby Isaac 562ad540459SPierre Jolivet for (j = 0; j < Nc * dimC; j++) rD[i * Nc * dimC + j] = 0.; 563805e7170SToby Isaac for (k = 0; k < Nq; k++) { 564805e7170SToby Isaac vol += geom->detJ[i * Nq + k] * weights[k]; 565805e7170SToby Isaac for (j = 0; j < Nc; j++) { 566805e7170SToby Isaac PetscReal pD[3] = {0., 0., 0.}; 567805e7170SToby Isaac 568805e7170SToby Isaac for (l = 0; l < dimC; l++) { 569ad540459SPierre 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]; 570805e7170SToby Isaac } 571ad540459SPierre Jolivet for (l = 0; l < dimC; l++) rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l]; 572805e7170SToby Isaac } 573805e7170SToby Isaac } 574805e7170SToby Isaac for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol; 575805e7170SToby Isaac } 576805e7170SToby Isaac } 577805e7170SToby Isaac } 578805e7170SToby Isaac if (H) { 579805e7170SToby Isaac PetscInt i, j, k, l, m, q, r; 580805e7170SToby Isaac 581805e7170SToby Isaac if (type == PETSC_SCALAR) { 582805e7170SToby Isaac PetscScalar *sH = (PetscScalar *)H; 583805e7170SToby Isaac PetscScalar *sqH = (PetscScalar *)qH; 584805e7170SToby Isaac 585805e7170SToby Isaac for (i = 0; i < numPoints; i++) { 586805e7170SToby Isaac PetscReal vol = 0.; 587805e7170SToby Isaac 588ad540459SPierre Jolivet for (j = 0; j < Nc * dimC * dimC; j++) sH[i * Nc * dimC * dimC + j] = 0.; 589805e7170SToby Isaac for (k = 0; k < Nq; k++) { 5904d1a973fSToby Isaac const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC]; 591805e7170SToby Isaac 592805e7170SToby Isaac vol += geom->detJ[i * Nq + k] * weights[k]; 593805e7170SToby Isaac for (j = 0; j < Nc; j++) { 5949371c9d4SSatish Balay PetscScalar pH[3][3] = { 5959371c9d4SSatish Balay {0., 0., 0.}, 5969371c9d4SSatish Balay {0., 0., 0.}, 5979371c9d4SSatish Balay {0., 0., 0.} 5989371c9d4SSatish Balay }; 599805e7170SToby Isaac const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC]; 600805e7170SToby Isaac 601805e7170SToby Isaac for (l = 0; l < dimC; l++) { 602805e7170SToby Isaac for (m = 0; m < dimC; m++) { 603805e7170SToby Isaac for (q = 0; q < dim; q++) { 604ad540459SPierre Jolivet for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r]; 605805e7170SToby Isaac } 606805e7170SToby Isaac } 607805e7170SToby Isaac } 608805e7170SToby Isaac for (l = 0; l < dimC; l++) { 609ad540459SPierre 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]; 610805e7170SToby Isaac } 611805e7170SToby Isaac } 612805e7170SToby Isaac } 613805e7170SToby Isaac for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol; 614805e7170SToby Isaac } 615805e7170SToby Isaac } else { 616805e7170SToby Isaac PetscReal *rH = (PetscReal *)H; 617805e7170SToby Isaac PetscReal *rqH = (PetscReal *)qH; 618805e7170SToby Isaac 619805e7170SToby Isaac for (i = 0; i < numPoints; i++) { 620805e7170SToby Isaac PetscReal vol = 0.; 621805e7170SToby Isaac 622ad540459SPierre Jolivet for (j = 0; j < Nc * dimC * dimC; j++) rH[i * Nc * dimC * dimC + j] = 0.; 623805e7170SToby Isaac for (k = 0; k < Nq; k++) { 624805e7170SToby Isaac const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC]; 625805e7170SToby Isaac 626805e7170SToby Isaac vol += geom->detJ[i * Nq + k] * weights[k]; 627805e7170SToby Isaac for (j = 0; j < Nc; j++) { 6289371c9d4SSatish Balay PetscReal pH[3][3] = { 6299371c9d4SSatish Balay {0., 0., 0.}, 6309371c9d4SSatish Balay {0., 0., 0.}, 6319371c9d4SSatish Balay {0., 0., 0.} 6329371c9d4SSatish Balay }; 633805e7170SToby Isaac const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC]; 634805e7170SToby Isaac 635805e7170SToby Isaac for (l = 0; l < dimC; l++) { 636805e7170SToby Isaac for (m = 0; m < dimC; m++) { 637805e7170SToby Isaac for (q = 0; q < dim; q++) { 638ad540459SPierre Jolivet for (r = 0; r < dim; r++) pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r]; 639805e7170SToby Isaac } 640805e7170SToby Isaac } 641805e7170SToby Isaac } 642805e7170SToby Isaac for (l = 0; l < dimC; l++) { 643ad540459SPierre 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]; 644805e7170SToby Isaac } 645805e7170SToby Isaac } 646805e7170SToby Isaac } 647805e7170SToby Isaac for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol; 648805e7170SToby Isaac } 649805e7170SToby Isaac } 650805e7170SToby Isaac } 6519566063dSJacob Faibussowitsch if (B) PetscCall(DMRestoreWorkArray(field->dm, N, mpitype, &qB)); 6529566063dSJacob Faibussowitsch if (D) PetscCall(DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD)); 6539566063dSJacob Faibussowitsch if (H) PetscCall(DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH)); 6549566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&geom)); 6559566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 657805e7170SToby Isaac } 658805e7170SToby Isaac 659d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree) 660d71ae5a4SJacob Faibussowitsch { 66151a454edSToby Isaac DMField_DS *dsfield; 66251a454edSToby Isaac PetscObject disc; 663741db25dSToby Isaac PetscInt h, imin, imax; 66451a454edSToby Isaac PetscClassId id; 66551a454edSToby Isaac 66651a454edSToby Isaac PetscFunctionBegin; 66751a454edSToby Isaac dsfield = (DMField_DS *)field->data; 6689566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(pointIS, &imin, &imax)); 669741db25dSToby Isaac if (imin >= imax) { 67051aa0bf7SToby Isaac h = 0; 67151aa0bf7SToby Isaac } else { 672f99c8401SToby Isaac for (h = 0; h < dsfield->height; h++) { 67351a454edSToby Isaac PetscInt hEnd; 67451a454edSToby Isaac 6759566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(field->dm, h, NULL, &hEnd)); 67651a454edSToby Isaac if (imin < hEnd) break; 67751a454edSToby Isaac } 67851aa0bf7SToby Isaac } 6796858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 68151a454edSToby Isaac if (id == PETSCFE_CLASSID) { 68251a454edSToby Isaac PetscFE fe = (PetscFE)disc; 68351a454edSToby Isaac PetscSpace sp; 68451a454edSToby Isaac 6859566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp)); 6869566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, minDegree, maxDegree)); 68751a454edSToby Isaac } 6883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68951a454edSToby Isaac } 69051a454edSToby Isaac 691d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldGetFVQuadrature_Internal(DMField field, IS pointIS, PetscQuadrature *quad) 692d71ae5a4SJacob Faibussowitsch { 693a2aba86cSMatthew G. Knepley DM dm = field->dm; 694a2aba86cSMatthew G. Knepley const PetscInt *points; 695a2aba86cSMatthew G. Knepley DMPolytopeType ct; 696a2aba86cSMatthew G. Knepley PetscInt dim, n; 697a2aba86cSMatthew G. Knepley PetscBool isplex; 698a2aba86cSMatthew G. Knepley 699a2aba86cSMatthew G. Knepley PetscFunctionBegin; 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex)); 7019566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &n)); 702a2aba86cSMatthew G. Knepley if (isplex && n) { 7039566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7049566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointIS, &points)); 7059566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, points[0], &ct)); 706a2aba86cSMatthew G. Knepley switch (ct) { 707a2aba86cSMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 708d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 709d71ae5a4SJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 1, -1.0, 1.0, quad)); 710d71ae5a4SJacob Faibussowitsch break; 711d71ae5a4SJacob Faibussowitsch default: 712d71ae5a4SJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad)); 713a2aba86cSMatthew G. Knepley } 7149566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointIS, &points)); 7151baa6e33SBarry Smith } else PetscCall(DMFieldCreateDefaultQuadrature(field, pointIS, quad)); 7163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 717a2aba86cSMatthew G. Knepley } 718a2aba86cSMatthew G. Knepley 719d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) 720d71ae5a4SJacob Faibussowitsch { 7211fdf1f07SMatthew G. Knepley PetscInt h, dim, imax, imin, cellHeight; 72251a454edSToby Isaac DM dm; 72351a454edSToby Isaac DMField_DS *dsfield; 72451a454edSToby Isaac PetscObject disc; 72551a454edSToby Isaac PetscFE fe; 72651a454edSToby Isaac PetscClassId id; 72751a454edSToby Isaac 72851a454edSToby Isaac PetscFunctionBegin; 72951a454edSToby Isaac dm = field->dm; 73051a454edSToby Isaac dsfield = (DMField_DS *)field->data; 7319566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(pointIS, &imin, &imax)); 7329566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 73351a454edSToby Isaac for (h = 0; h <= dim; h++) { 73451a454edSToby Isaac PetscInt hStart, hEnd; 73551a454edSToby Isaac 7369566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd)); 7376c041b70SSatish Balay if (imax >= hStart && imin < hEnd) break; 73851a454edSToby Isaac } 7399566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 7401fdf1f07SMatthew G. Knepley h -= cellHeight; 74151a454edSToby Isaac *quad = NULL; 742f99c8401SToby Isaac if (h < dsfield->height) { 7436858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc)); 7449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 7453ba16761SJacob Faibussowitsch if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS); 74651a454edSToby Isaac fe = (PetscFE)disc; 7479566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, quad)); 7489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)*quad)); 74951a454edSToby Isaac } 7503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75151a454edSToby Isaac } 75251a454edSToby Isaac 753989fa639SBrad Aagaard static PetscErrorCode DMFieldCreateDefaultFaceQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad) 754989fa639SBrad Aagaard { 755989fa639SBrad Aagaard PetscInt h, dim, imax, imin, cellHeight; 756989fa639SBrad Aagaard DM dm; 757989fa639SBrad Aagaard DMField_DS *dsfield; 758989fa639SBrad Aagaard PetscObject disc; 759989fa639SBrad Aagaard PetscFE fe; 760989fa639SBrad Aagaard PetscClassId id; 761989fa639SBrad Aagaard 762989fa639SBrad Aagaard PetscFunctionBegin; 763989fa639SBrad Aagaard dm = field->dm; 764989fa639SBrad Aagaard dsfield = (DMField_DS *)field->data; 765989fa639SBrad Aagaard PetscCall(ISGetMinMax(pointIS, &imin, &imax)); 766989fa639SBrad Aagaard PetscCall(DMGetDimension(dm, &dim)); 767989fa639SBrad Aagaard for (h = 0; h <= dim; h++) { 768989fa639SBrad Aagaard PetscInt hStart, hEnd; 769989fa639SBrad Aagaard 770989fa639SBrad Aagaard PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd)); 771989fa639SBrad Aagaard if (imax >= hStart && imin < hEnd) break; 772989fa639SBrad Aagaard } 773989fa639SBrad Aagaard PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 774989fa639SBrad Aagaard h -= cellHeight; 775989fa639SBrad Aagaard *quad = NULL; 776989fa639SBrad Aagaard if (h < dsfield->height) { 777989fa639SBrad Aagaard PetscQuadrature fq; 778989fa639SBrad Aagaard 779989fa639SBrad Aagaard PetscCall(DMFieldDSGetHeightDisc(field, h, dsfield->disc, &disc)); 780989fa639SBrad Aagaard PetscCall(PetscObjectGetClassId(disc, &id)); 781989fa639SBrad Aagaard if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS); 782989fa639SBrad Aagaard fe = (PetscFE)disc; 783989fa639SBrad Aagaard PetscCall(PetscFEGetFaceQuadrature(fe, &fq)); 784989fa639SBrad Aagaard PetscCall(PetscFEExpandFaceQuadrature(fe, fq, quad)); 785989fa639SBrad Aagaard } 786989fa639SBrad Aagaard PetscFunctionReturn(PETSC_SUCCESS); 787989fa639SBrad Aagaard } 788989fa639SBrad Aagaard 789d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom) 790d71ae5a4SJacob Faibussowitsch { 7910145028aSToby Isaac const PetscInt *points; 7920145028aSToby Isaac PetscInt p, dim, dE, numFaces, Nq; 793b7260050SToby Isaac PetscInt maxDegree; 7940145028aSToby Isaac DMLabel depthLabel; 7950145028aSToby Isaac IS cellIS; 7960145028aSToby Isaac DM dm = field->dm; 7970145028aSToby Isaac 7980145028aSToby Isaac PetscFunctionBegin; 7990145028aSToby Isaac dim = geom->dim; 8000145028aSToby Isaac dE = geom->dimEmbed; 8019566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depthLabel)); 8029566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS)); 8039566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(field, cellIS, NULL, &maxDegree)); 8049566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointIS, &points)); 8050145028aSToby Isaac numFaces = geom->numCells; 8060145028aSToby Isaac Nq = geom->numPoints; 807f15274beSMatthew Knepley /* First, set local faces and flip normals so that they are outward for the first supporting cell */ 808f15274beSMatthew Knepley for (p = 0; p < numFaces; p++) { 809f15274beSMatthew Knepley PetscInt point = points[p]; 810f15274beSMatthew Knepley PetscInt suppSize, s, coneSize, c, numChildren; 8110e18dc48SMatthew G. Knepley const PetscInt *supp; 812f15274beSMatthew Knepley 8139566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL)); 81428b400f6SJacob Faibussowitsch PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 8159566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, point, &suppSize)); 81663a3b9bcSJacob Faibussowitsch PetscCheck(suppSize <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, suppSize); 817f15274beSMatthew Knepley if (!suppSize) continue; 8189566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, point, &supp)); 819f15274beSMatthew Knepley for (s = 0; s < suppSize; ++s) { 8200e18dc48SMatthew G. Knepley const PetscInt *cone, *ornt; 8210e18dc48SMatthew G. Knepley 8229566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, supp[s], &coneSize)); 8230e18dc48SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, supp[s], &cone, &ornt)); 8249371c9d4SSatish Balay for (c = 0; c < coneSize; ++c) 8259371c9d4SSatish Balay if (cone[c] == point) break; 8260e18dc48SMatthew G. Knepley geom->face[p][s * 2 + 0] = c; 8270e18dc48SMatthew G. Knepley geom->face[p][s * 2 + 1] = ornt[c]; 82884cbd072SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, supp[s], &cone, &ornt)); 82984cbd072SMatthew 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]); 830f15274beSMatthew Knepley } 8310e18dc48SMatthew G. Knepley if (geom->face[p][1] < 0) { 832f15274beSMatthew Knepley PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d; 833f15274beSMatthew Knepley 8349371c9d4SSatish Balay for (q = 0; q < Np; ++q) 8359371c9d4SSatish Balay for (d = 0; d < dE; ++d) geom->n[(p * Np + q) * dE + d] = -geom->n[(p * Np + q) * dE + d]; 836f15274beSMatthew Knepley } 837f15274beSMatthew Knepley } 838b7260050SToby Isaac if (maxDegree <= 1) { 839313f16f2SMatthew G. Knepley PetscQuadrature cellQuad = NULL; 8400145028aSToby Isaac PetscInt numCells, offset, *cells; 8410145028aSToby Isaac PetscFEGeom *cellGeom; 8420145028aSToby Isaac IS suppIS; 8430145028aSToby Isaac 844313f16f2SMatthew G. Knepley if (quad) { 845313f16f2SMatthew G. Knepley DM dm; 846313f16f2SMatthew G. Knepley PetscReal *points, *weights; 847313f16f2SMatthew G. Knepley PetscInt tdim, Nc, Np; 848313f16f2SMatthew G. Knepley 849313f16f2SMatthew G. Knepley PetscCall(DMFieldGetDM(field, &dm)); 850313f16f2SMatthew G. Knepley PetscCall(DMGetDimension(dm, &tdim)); 851313f16f2SMatthew G. Knepley if (tdim > dim) { 852313f16f2SMatthew G. Knepley // Make a compatible cell quadrature (points don't matter since its affine) 853313f16f2SMatthew G. Knepley PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad)); 854313f16f2SMatthew G. Knepley PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Np, NULL, NULL)); 855313f16f2SMatthew G. Knepley PetscCall(PetscCalloc1((dim + 1) * Np, &points)); 856313f16f2SMatthew G. Knepley PetscCall(PetscCalloc1(Nc * Np, &weights)); 857313f16f2SMatthew G. Knepley PetscCall(PetscQuadratureSetData(cellQuad, dim + 1, Nc, Np, points, weights)); 858313f16f2SMatthew G. Knepley } else { 859313f16f2SMatthew G. Knepley // TODO J will be wrong here, but other things need to be fixed 860313f16f2SMatthew G. Knepley // This path comes from calling DMProjectBdFieldLabelLocal() in Plex ex5 861313f16f2SMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)quad)); 862313f16f2SMatthew G. Knepley cellQuad = quad; 863313f16f2SMatthew G. Knepley } 864313f16f2SMatthew G. Knepley } 8650145028aSToby Isaac for (p = 0, numCells = 0; p < numFaces; p++) { 8660145028aSToby Isaac PetscInt point = points[p]; 8670145028aSToby Isaac PetscInt numSupp, numChildren; 8680145028aSToby Isaac 8699566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL)); 87028b400f6SJacob Faibussowitsch PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 8719566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, point, &numSupp)); 87263a3b9bcSJacob Faibussowitsch PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp); 8730145028aSToby Isaac numCells += numSupp; 8740145028aSToby Isaac } 8759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, &cells)); 8760145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 8770145028aSToby Isaac PetscInt point = points[p]; 8780145028aSToby Isaac PetscInt numSupp, s; 8790145028aSToby Isaac const PetscInt *supp; 8800145028aSToby Isaac 8819566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, point, &numSupp)); 8829566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, point, &supp)); 883ad540459SPierre Jolivet for (s = 0; s < numSupp; s++, offset++) cells[offset] = supp[s]; 8840145028aSToby Isaac } 8859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS)); 886ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FEGEOM_BASIC, &cellGeom)); 8870145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 8880145028aSToby Isaac PetscInt point = points[p]; 8890145028aSToby Isaac PetscInt numSupp, s, q; 8900145028aSToby Isaac const PetscInt *supp; 8910145028aSToby Isaac 8929566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, point, &numSupp)); 8939566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, point, &supp)); 8940145028aSToby Isaac for (s = 0; s < numSupp; s++, offset++) { 8950145028aSToby Isaac for (q = 0; q < Nq * dE * dE; q++) { 89641418987SMatthew G. Knepley geom->suppJ[s][p * Nq * dE * dE + q] = cellGeom->J[offset * Nq * dE * dE + q]; 8970145028aSToby Isaac geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 8980145028aSToby Isaac } 89941418987SMatthew G. Knepley for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q]; 9000145028aSToby Isaac } 9010145028aSToby Isaac } 9029566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&cellGeom)); 903313f16f2SMatthew G. Knepley PetscCall(PetscQuadratureDestroy(&cellQuad)); 9049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&suppIS)); 9059566063dSJacob Faibussowitsch PetscCall(PetscFree(cells)); 9060145028aSToby Isaac } else { 9076858538eSMatthew G. Knepley DMField_DS *dsfield = (DMField_DS *)field->data; 9080145028aSToby Isaac PetscObject faceDisc, cellDisc; 9090145028aSToby Isaac PetscClassId faceId, cellId; 9100145028aSToby Isaac PetscDualSpace dsp; 9110145028aSToby Isaac DM K; 912ba2698f1SMatthew G. Knepley DMPolytopeType ct; 9130145028aSToby Isaac PetscInt (*co)[2][3]; 9140145028aSToby Isaac PetscInt coneSize; 9150145028aSToby Isaac PetscInt **counts; 9160145028aSToby Isaac PetscInt f, i, o, q, s; 917e59eb9b3SMatthew G. Knepley PetscBool found = PETSC_FALSE; 9180145028aSToby Isaac const PetscInt *coneK; 919ba2698f1SMatthew G. Knepley PetscInt eStart, minOrient, maxOrient, numOrient; 9200145028aSToby Isaac PetscInt *orients; 9210145028aSToby Isaac PetscReal **orientPoints; 9220145028aSToby Isaac PetscReal *cellPoints; 9230145028aSToby Isaac PetscReal *dummyWeights; 9240145028aSToby Isaac PetscQuadrature cellQuad = NULL; 9250145028aSToby Isaac 9266858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, 1, dsfield->disc, &faceDisc)); 9276858538eSMatthew G. Knepley PetscCall(DMFieldDSGetHeightDisc(field, 0, dsfield->disc, &cellDisc)); 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(faceDisc, &faceId)); 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(cellDisc, &cellId)); 9301dca8a05SBarry Smith PetscCheck(faceId == PETSCFE_CLASSID && cellId == PETSCFE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported"); 9319566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace((PetscFE)cellDisc, &dsp)); 9329566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(dsp, &K)); 9339566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(K, 1, &eStart, NULL)); 9349566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(K, eStart, &ct)); 9359566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(K, 0, &coneSize)); 9369566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(K, 0, &coneK)); 9379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numFaces, &co, coneSize, &counts)); 9389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dE * Nq, &cellPoints)); 9399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nq, &dummyWeights)); 9409566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad)); 9419566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights)); 9421690c2aeSBarry Smith minOrient = PETSC_INT_MAX; 9431690c2aeSBarry Smith maxOrient = PETSC_INT_MIN; 9440145028aSToby Isaac for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */ 9450145028aSToby Isaac PetscInt point = points[p]; 9460145028aSToby Isaac PetscInt numSupp, numChildren; 9470145028aSToby Isaac const PetscInt *supp; 9480145028aSToby Isaac 9499566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, point, &numChildren, NULL)); 95028b400f6SJacob Faibussowitsch PetscCheck(!numChildren, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children"); 9519566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, point, &numSupp)); 95263a3b9bcSJacob Faibussowitsch PetscCheck(numSupp <= 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %" PetscInt_FMT " has %" PetscInt_FMT " support, expected at most 2", point, numSupp); 9539566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, point, &supp)); 9540145028aSToby Isaac for (s = 0; s < numSupp; s++) { 9550145028aSToby Isaac PetscInt cell = supp[s]; 9560145028aSToby Isaac PetscInt numCone; 9570145028aSToby Isaac const PetscInt *cone, *orient; 9580145028aSToby Isaac 9599566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &numCone)); 960e59eb9b3SMatthew G. Knepley // When we extract submeshes, we hang cells from the side that are not fully realized. We ignore these 961e59eb9b3SMatthew G. Knepley if (numCone == 1) { 962e59eb9b3SMatthew G. Knepley co[p][s][0] = -1; 963e59eb9b3SMatthew G. Knepley co[p][s][1] = -1; 964e59eb9b3SMatthew G. Knepley co[p][s][2] = -1; 965e59eb9b3SMatthew G. Knepley continue; 966e59eb9b3SMatthew G. Knepley } 96708401ef6SPierre Jolivet PetscCheck(numCone == coneSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element"); 9689566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 9699566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cell, &orient)); 9700145028aSToby Isaac for (f = 0; f < coneSize; f++) { 9710145028aSToby Isaac if (cone[f] == point) break; 9720145028aSToby Isaac } 9730145028aSToby Isaac co[p][s][0] = f; 9740145028aSToby Isaac co[p][s][1] = orient[f]; 9750145028aSToby Isaac co[p][s][2] = cell; 9760145028aSToby Isaac minOrient = PetscMin(minOrient, orient[f]); 977cd93b0e1SMatthew G. Knepley maxOrient = PetscMax(maxOrient, orient[f]); 978e59eb9b3SMatthew G. Knepley found = PETSC_TRUE; 9790145028aSToby Isaac } 9800145028aSToby Isaac for (; s < 2; s++) { 9810145028aSToby Isaac co[p][s][0] = -1; 9820145028aSToby Isaac co[p][s][1] = -1; 9830145028aSToby Isaac co[p][s][2] = -1; 9840145028aSToby Isaac } 9850145028aSToby Isaac } 986e59eb9b3SMatthew G. Knepley numOrient = found ? maxOrient + 1 - minOrient : 0; 9879566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(K, 0, &coneK)); 9880145028aSToby Isaac /* count all (face,orientation) doubles that appear */ 9899566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numOrient, &orients, numOrient, &orientPoints)); 9909566063dSJacob Faibussowitsch for (f = 0; f < coneSize; f++) PetscCall(PetscCalloc1(numOrient + 1, &counts[f])); 9910145028aSToby Isaac for (p = 0; p < numFaces; p++) { 9920145028aSToby Isaac for (s = 0; s < 2; s++) { 9930145028aSToby Isaac if (co[p][s][0] >= 0) { 9940145028aSToby Isaac counts[co[p][s][0]][co[p][s][1] - minOrient]++; 9950145028aSToby Isaac orients[co[p][s][1] - minOrient]++; 9960145028aSToby Isaac } 9970145028aSToby Isaac } 9980145028aSToby Isaac } 9990145028aSToby Isaac for (o = 0; o < numOrient; o++) { 10000145028aSToby Isaac if (orients[o]) { 10010145028aSToby Isaac PetscInt orient = o + minOrient; 10020145028aSToby Isaac PetscInt q; 10030145028aSToby Isaac 10049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nq * dim, &orientPoints[o])); 10050145028aSToby Isaac /* rotate the quadrature points appropriately */ 1006ba2698f1SMatthew G. Knepley switch (ct) { 1007d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 1008d71ae5a4SJacob Faibussowitsch break; 1009ba2698f1SMatthew G. Knepley case DM_POLYTOPE_SEGMENT: 10100145028aSToby Isaac if (orient == -2 || orient == 1) { 1011ad540459SPierre Jolivet for (q = 0; q < Nq; q++) orientPoints[o][q] = -geom->xi[q]; 10120145028aSToby Isaac } else { 1013ad540459SPierre Jolivet for (q = 0; q < Nq; q++) orientPoints[o][q] = geom->xi[q]; 10140145028aSToby Isaac } 10150145028aSToby Isaac break; 1016ba2698f1SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 10170145028aSToby Isaac for (q = 0; q < Nq; q++) { 10180145028aSToby Isaac PetscReal lambda[3]; 10190145028aSToby Isaac PetscReal lambdao[3]; 10200145028aSToby Isaac 10210145028aSToby Isaac /* convert to barycentric */ 10220145028aSToby Isaac lambda[0] = -(geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.; 10230145028aSToby Isaac lambda[1] = (geom->xi[2 * q] + 1.) / 2.; 10240145028aSToby Isaac lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.; 10250145028aSToby Isaac if (orient >= 0) { 1026ad540459SPierre Jolivet for (i = 0; i < 3; i++) lambdao[i] = lambda[(orient + i) % 3]; 10270145028aSToby Isaac } else { 1028ad540459SPierre Jolivet for (i = 0; i < 3; i++) lambdao[i] = lambda[(-(orient + i) + 3) % 3]; 10290145028aSToby Isaac } 10300145028aSToby Isaac /* convert to coordinates */ 10310145028aSToby Isaac orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1]; 10320145028aSToby Isaac orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2]; 10330145028aSToby Isaac } 10340145028aSToby Isaac break; 1035ba2698f1SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 10360145028aSToby Isaac for (q = 0; q < Nq; q++) { 10370145028aSToby Isaac PetscReal xi[2], xio[2]; 10380145028aSToby Isaac PetscInt oabs = (orient >= 0) ? orient : -(orient + 1); 10390145028aSToby Isaac 10400145028aSToby Isaac xi[0] = geom->xi[2 * q]; 10410145028aSToby Isaac xi[1] = geom->xi[2 * q + 1]; 10420145028aSToby Isaac switch (oabs) { 10430145028aSToby Isaac case 1: 10440145028aSToby Isaac xio[0] = xi[1]; 10450145028aSToby Isaac xio[1] = -xi[0]; 10460145028aSToby Isaac break; 1047f4d061e9SPierre Jolivet case 2: 1048f4d061e9SPierre Jolivet xio[0] = -xi[0]; 1049f4d061e9SPierre Jolivet xio[1] = -xi[1]; 1050f4d061e9SPierre Jolivet break; 1051f4d061e9SPierre Jolivet case 3: 1052f4d061e9SPierre Jolivet xio[0] = -xi[1]; 1053f4d061e9SPierre Jolivet xio[1] = xi[0]; 1054f4d061e9SPierre Jolivet break; 1055aa1371cbSToby Isaac case 0: 1056aa1371cbSToby Isaac default: 1057aa1371cbSToby Isaac xio[0] = xi[0]; 1058aa1371cbSToby Isaac xio[1] = xi[1]; 1059aa1371cbSToby Isaac break; 10600145028aSToby Isaac } 1061ad540459SPierre Jolivet if (orient < 0) xio[0] = -xio[0]; 10620145028aSToby Isaac orientPoints[o][2 * q + 0] = xio[0]; 10630145028aSToby Isaac orientPoints[o][2 * q + 1] = xio[1]; 10640145028aSToby Isaac } 10650145028aSToby Isaac break; 1066d71ae5a4SJacob Faibussowitsch default: 1067d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s not yet supported", DMPolytopeTypes[ct]); 10680145028aSToby Isaac } 10690145028aSToby Isaac } 10700145028aSToby Isaac } 10710145028aSToby Isaac for (f = 0; f < coneSize; f++) { 10720145028aSToby Isaac PetscInt face = coneK[f]; 10734d1a973fSToby Isaac PetscReal v0[3]; 10742d89661fSMatthew G. Knepley PetscReal J[9], detJ; 10750145028aSToby Isaac PetscInt numCells, offset; 10760145028aSToby Isaac PetscInt *cells; 10770145028aSToby Isaac IS suppIS; 10780145028aSToby Isaac 10799566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ)); 10800145028aSToby Isaac for (o = 0; o <= numOrient; o++) { 10810145028aSToby Isaac PetscFEGeom *cellGeom; 10820145028aSToby Isaac 10830145028aSToby Isaac if (!counts[f][o]) continue; 10840145028aSToby Isaac /* If this (face,orientation) double appears, 10850145028aSToby Isaac * convert the face quadrature points into volume quadrature points */ 10860145028aSToby Isaac for (q = 0; q < Nq; q++) { 10870145028aSToby Isaac PetscReal xi0[3] = {-1., -1., -1.}; 10880145028aSToby Isaac 10890145028aSToby Isaac CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]); 10900145028aSToby Isaac } 10910145028aSToby Isaac for (p = 0, numCells = 0; p < numFaces; p++) { 10920145028aSToby Isaac for (s = 0; s < 2; s++) { 10930145028aSToby Isaac if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++; 10940145028aSToby Isaac } 10950145028aSToby Isaac } 10969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, &cells)); 10970145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 10980145028aSToby Isaac for (s = 0; s < 2; s++) { 1099ad540459SPierre Jolivet if (co[p][s][0] == f && co[p][s][1] == o + minOrient) cells[offset++] = co[p][s][2]; 11000145028aSToby Isaac } 11010145028aSToby Isaac } 11029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numCells, cells, PETSC_USE_POINTER, &suppIS)); 1103ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(field, suppIS, cellQuad, PETSC_FEGEOM_BASIC, &cellGeom)); 11040145028aSToby Isaac for (p = 0, offset = 0; p < numFaces; p++) { 11050145028aSToby Isaac for (s = 0; s < 2; s++) { 11060145028aSToby Isaac if (co[p][s][0] == f && co[p][s][1] == o + minOrient) { 11070145028aSToby Isaac for (q = 0; q < Nq * dE * dE; q++) { 110841418987SMatthew G. Knepley geom->suppJ[s][p * Nq * dE * dE + q] = cellGeom->J[offset * Nq * dE * dE + q]; 11090145028aSToby Isaac geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q]; 11100145028aSToby Isaac } 111141418987SMatthew G. Knepley for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q]; 11120145028aSToby Isaac offset++; 11130145028aSToby Isaac } 11140145028aSToby Isaac } 11150145028aSToby Isaac } 11169566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&cellGeom)); 11179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&suppIS)); 11189566063dSJacob Faibussowitsch PetscCall(PetscFree(cells)); 11190145028aSToby Isaac } 11200145028aSToby Isaac } 11210145028aSToby Isaac for (o = 0; o < numOrient; o++) { 112248a46eb9SPierre Jolivet if (orients[o]) PetscCall(PetscFree(orientPoints[o])); 11230145028aSToby Isaac } 11249566063dSJacob Faibussowitsch PetscCall(PetscFree2(orients, orientPoints)); 11259566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&cellQuad)); 11269566063dSJacob Faibussowitsch for (f = 0; f < coneSize; f++) PetscCall(PetscFree(counts[f])); 11279566063dSJacob Faibussowitsch PetscCall(PetscFree2(co, counts)); 11280145028aSToby Isaac } 11299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointIS, &points)); 11309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 11313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11320145028aSToby Isaac } 11330145028aSToby Isaac 1134d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMFieldInitialize_DS(DMField field) 1135d71ae5a4SJacob Faibussowitsch { 113651a454edSToby Isaac PetscFunctionBegin; 113751a454edSToby Isaac field->ops->destroy = DMFieldDestroy_DS; 113851a454edSToby Isaac field->ops->evaluate = DMFieldEvaluate_DS; 113951a454edSToby Isaac field->ops->evaluateFE = DMFieldEvaluateFE_DS; 1140805e7170SToby Isaac field->ops->evaluateFV = DMFieldEvaluateFV_DS; 1141b7260050SToby Isaac field->ops->getDegree = DMFieldGetDegree_DS; 114251a454edSToby Isaac field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS; 1143989fa639SBrad Aagaard field->ops->createDefaultFaceQuadrature = DMFieldCreateDefaultFaceQuadrature_DS; 114451a454edSToby Isaac field->ops->view = DMFieldView_DS; 11450145028aSToby Isaac field->ops->computeFaceData = DMFieldComputeFaceData_DS; 11463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114751a454edSToby Isaac } 114851a454edSToby Isaac 1149d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field) 1150d71ae5a4SJacob Faibussowitsch { 115151a454edSToby Isaac DMField_DS *dsfield; 115251a454edSToby Isaac 115351a454edSToby Isaac PetscFunctionBegin; 11544dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dsfield)); 115551a454edSToby Isaac field->data = dsfield; 11569566063dSJacob Faibussowitsch PetscCall(DMFieldInitialize_DS(field)); 11573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115851a454edSToby Isaac } 115951a454edSToby Isaac 1160d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDSWithDG(DM dm, DM dmDG, PetscInt fieldNum, Vec vec, Vec vecDG, DMField *field) 1161d71ae5a4SJacob Faibussowitsch { 116251a454edSToby Isaac DMField b; 116351a454edSToby Isaac DMField_DS *dsfield; 11646858538eSMatthew G. Knepley PetscObject disc = NULL, discDG = NULL; 11656858538eSMatthew G. Knepley PetscSection section; 116651a454edSToby Isaac PetscBool isContainer = PETSC_FALSE; 116751a454edSToby Isaac PetscClassId id = -1; 116851a454edSToby Isaac PetscInt numComponents = -1, dsNumFields; 116951a454edSToby Isaac 117051a454edSToby Isaac PetscFunctionBegin; 11719566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 11729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(section, fieldNum, &numComponents)); 11739566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &dsNumFields)); 11749566063dSJacob Faibussowitsch if (dsNumFields) PetscCall(DMGetField(dm, fieldNum, NULL, &disc)); 11756858538eSMatthew G. Knepley if (dsNumFields && dmDG) { 11766858538eSMatthew G. Knepley PetscCall(DMGetField(dmDG, fieldNum, NULL, &discDG)); 11776858538eSMatthew G. Knepley PetscCall(PetscObjectReference(discDG)); 11786858538eSMatthew G. Knepley } 117951a454edSToby Isaac if (disc) { 11809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 118151a454edSToby Isaac isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE; 118251a454edSToby Isaac } 118351a454edSToby Isaac if (!disc || isContainer) { 118451a454edSToby Isaac MPI_Comm comm = PetscObjectComm((PetscObject)dm); 118551a454edSToby Isaac PetscFE fe; 11862df84da0SMatthew G. Knepley DMPolytopeType ct, locct = DM_POLYTOPE_UNKNOWN; 11872df84da0SMatthew G. Knepley PetscInt dim, cStart, cEnd, cellHeight; 118851a454edSToby Isaac 11899566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 11909566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 11919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 11929566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &locct)); 11936a210b70SBarry Smith PetscCallMPI(MPIU_Allreduce(&locct, &ct, 1, MPI_INT, MPI_MIN, comm)); 11949566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, numComponents, ct, 1, PETSC_DETERMINE, &fe)); 11959566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(fe, NULL, "-field_fe_view")); 119651a454edSToby Isaac disc = (PetscObject)fe; 11971baa6e33SBarry Smith } else PetscCall(PetscObjectReference(disc)); 11989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 11996858538eSMatthew G. Knepley if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)disc, &numComponents)); 12006858538eSMatthew G. Knepley else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine number of discretization components"); 12019566063dSJacob Faibussowitsch PetscCall(DMFieldCreate(dm, numComponents, DMFIELD_VERTEX, &b)); 12029566063dSJacob Faibussowitsch PetscCall(DMFieldSetType(b, DMFIELDDS)); 120351a454edSToby Isaac dsfield = (DMField_DS *)b->data; 120451a454edSToby Isaac dsfield->fieldNum = fieldNum; 12059566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dsfield->height)); 1206f99c8401SToby Isaac dsfield->height++; 12079566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dsfield->height, &dsfield->disc)); 120851a454edSToby Isaac dsfield->disc[0] = disc; 12099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)vec)); 121051a454edSToby Isaac dsfield->vec = vec; 12116858538eSMatthew G. Knepley if (dmDG) { 12126858538eSMatthew G. Knepley dsfield->dmDG = dmDG; 12136858538eSMatthew G. Knepley PetscCall(PetscCalloc1(dsfield->height, &dsfield->discDG)); 12146858538eSMatthew G. Knepley dsfield->discDG[0] = discDG; 12156858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)vecDG)); 12166858538eSMatthew G. Knepley dsfield->vecDG = vecDG; 12176858538eSMatthew G. Knepley } 121851a454edSToby Isaac *field = b; 12193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12206858538eSMatthew G. Knepley } 122151a454edSToby Isaac 1222d71ae5a4SJacob Faibussowitsch PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec, DMField *field) 1223d71ae5a4SJacob Faibussowitsch { 12246858538eSMatthew G. Knepley PetscFunctionBegin; 12256858538eSMatthew G. Knepley PetscCall(DMFieldCreateDSWithDG(dm, NULL, fieldNum, vec, NULL, field)); 12263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122751a454edSToby Isaac } 1228