xref: /petsc/src/dm/impls/plex/plexproject.c (revision 69291d52fa0983ddcb3e907357ab97df51e3b8d8)
147923291SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
247923291SMatthew G. Knepley 
38c6c5593SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
48c6c5593SMatthew G. Knepley 
5496733e2SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscReal time, PetscFECellGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
68c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
78c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
847923291SMatthew G. Knepley {
98c6c5593SMatthew G. Knepley   PetscDS        prob;
109c3cf19fSMatthew G. Knepley   PetscInt       Nf, *Nc, f, totDim, spDim, d, v;
118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
128c6c5593SMatthew G. Knepley 
138c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
148c6c5593SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
158c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
16496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
179c3cf19fSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
188c6c5593SMatthew G. Knepley   /* Get values for closure */
198c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
208c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
218c6c5593SMatthew G. Knepley 
228c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
238c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
249c3cf19fSMatthew G. Knepley     for (d = 0; d < spDim; ++d, ++v) {
258c6c5593SMatthew G. Knepley       if (funcs[f]) {
268c6c5593SMatthew G. Knepley         if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, fegeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
278c6c5593SMatthew G. Knepley         else         {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
288c6c5593SMatthew G. Knepley       } else {
299c3cf19fSMatthew G. Knepley         values[v] = 0.0;
308c6c5593SMatthew G. Knepley       }
318c6c5593SMatthew G. Knepley     }
329c3cf19fSMatthew G. Knepley   }
338c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
348c6c5593SMatthew G. Knepley }
358c6c5593SMatthew G. Knepley 
361c531cf8SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, DM dmAux, PetscReal time, Vec localU, Vec localA, PetscFECellGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
370c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
388c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
398c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
408c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
41191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
428c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
438c6c5593SMatthew G. Knepley {
448c6c5593SMatthew G. Knepley   PetscDS            prob, probAux = NULL;
458c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
469c3cf19fSMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
478c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
488c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
49191494d9SMatthew G. Knepley   const PetscScalar *constants;
508c6c5593SMatthew G. Knepley   PetscReal         *x;
51496733e2SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
5216457704SMatthew G. Knepley   PetscInt           dim, dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, c, v, tp = 0;
538c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
548c6c5593SMatthew G. Knepley 
558c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
568c6c5593SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
578c6c5593SMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
588c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
59496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
60496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
618c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
628c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
639c3cf19fSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
64496733e2SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
65191494d9SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
668c6c5593SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
67496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
688c6c5593SMatthew G. Knepley   if (dmAux) {
691ba34526SMatthew G. Knepley     DMLabel  spmap;
701ba34526SMatthew G. Knepley     PetscInt subp;
711ba34526SMatthew G. Knepley 
721ba34526SMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
731ba34526SMatthew G. Knepley     if (spmap) {
741ba34526SMatthew G. Knepley       IS              subpointIS;
751ba34526SMatthew G. Knepley       const PetscInt *subpoints;
761ba34526SMatthew G. Knepley       PetscInt        numSubpoints;
771ba34526SMatthew G. Knepley 
781ba34526SMatthew G. Knepley       ierr = DMPlexCreateSubpointIS(dmAux, &subpointIS);CHKERRQ(ierr);
791ba34526SMatthew G. Knepley       ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr);
801ba34526SMatthew G. Knepley       ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr);
811ba34526SMatthew G. Knepley       ierr = PetscFindInt(p, numSubpoints, subpoints, &subp);CHKERRQ(ierr);
821ba34526SMatthew G. Knepley       if (subp < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d not found in submesh", p);
831ba34526SMatthew G. Knepley       ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);
841ba34526SMatthew G. Knepley       ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
851ba34526SMatthew G. Knepley     } else subp = p;
868c6c5593SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
871ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
888c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
89496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
90496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
918c6c5593SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
928c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
938c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
9411d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
95496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
961ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
978c6c5593SMatthew G. Knepley   }
988c6c5593SMatthew G. Knepley   /* Get values for closure */
998c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
1008c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1018c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
1029c3cf19fSMatthew G. Knepley     for (d = 0; d < spDim; ++d, ++v) {
1038c6c5593SMatthew G. Knepley       if (funcs[f]) {
1048c6c5593SMatthew G. Knepley         PetscQuadrature  quad;
1059c3cf19fSMatthew G. Knepley         const PetscReal *points, *weights;
1068c6c5593SMatthew G. Knepley         PetscInt         numPoints, q;
1078c6c5593SMatthew G. Knepley 
1088c6c5593SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
1099c3cf19fSMatthew G. Knepley         ierr = PetscQuadratureGetData(quad, NULL, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
1100c898c18SMatthew G. Knepley         for (q = 0; q < numPoints; ++q, ++tp) {
1118c6c5593SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, fegeom->v0, fegeom->J, &points[q*dim], x);
112496733e2SMatthew G. Knepley           EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t);
1131ba34526SMatthew G. Knepley           if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
114191494d9SMatthew G. Knepley           (*funcs[f])(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, x, numConstants, constants, bc);
1151c531cf8SMatthew G. Knepley           if (Ncc) {
1161c531cf8SMatthew G. Knepley             for (c = 0; c < Ncc; ++c) values[v] += bc[comps[c]]*weights[comps[c]];
1171c531cf8SMatthew G. Knepley           } else {
1189c3cf19fSMatthew G. Knepley             for (c = 0; c < Nc[f]; ++c) values[v] += bc[c]*weights[c];
1198c6c5593SMatthew G. Knepley           }
1201c531cf8SMatthew G. Knepley         }
1218c6c5593SMatthew G. Knepley       } else {
1229c3cf19fSMatthew G. Knepley         values[v] = 0.0;
1238c6c5593SMatthew G. Knepley       }
1248c6c5593SMatthew G. Knepley     }
1258c6c5593SMatthew G. Knepley   }
1268c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1278c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1288c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1298c6c5593SMatthew G. Knepley }
1308c6c5593SMatthew G. Knepley 
1310c898c18SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Private(DM dm, DM dmAux, PetscInt h, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
1321c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
1331c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1348c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1358c6c5593SMatthew G. Knepley {
1368c6c5593SMatthew G. Knepley   PetscFECellGeom fegeom;
1378c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1388c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1398c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1408c6c5593SMatthew G. Knepley 
1418c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1428c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1438c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1448c6c5593SMatthew G. Knepley   if (hasFE) {
14543048b12SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr);
1468c6c5593SMatthew G. Knepley     fegeom.dim      = dim - h;
1478c6c5593SMatthew G. Knepley     fegeom.dimEmbed = dimEmbed;
1488c6c5593SMatthew G. Knepley   }
1498c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1508c6c5593SMatthew G. Knepley   switch (type) {
1518c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1528c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
153496733e2SMatthew G. Knepley     ierr = DMProjectPoint_Func_Private(dm, time, &fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
1548c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1558c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
1561c531cf8SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, dmAux, time, localU, localA, &fegeom, sp, p, Ncc, comps,
1570c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1580c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1590c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1600c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
161191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1628c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1638c6c5593SMatthew G. Knepley   }
1648c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1658c6c5593SMatthew G. Knepley }
1668c6c5593SMatthew G. Knepley 
16746fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
1681c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
1698c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
1708c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
1718c6c5593SMatthew G. Knepley {
1728c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
1730c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
1748c6c5593SMatthew G. Knepley   Vec             localA = NULL;
17547923291SMatthew G. Knepley   PetscSection    section;
1768c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
1770c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
1788c6c5593SMatthew G. Knepley   PetscInt       *Nc;
1792af58022SMatthew G. Knepley   PetscInt        dim, dimEmbed, minHeight, maxHeight, h, Nf, NfAux = 0, f;
18088aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
18147923291SMatthew G. Knepley   PetscErrorCode  ierr;
18247923291SMatthew G. Knepley 
18347923291SMatthew G. Knepley   PetscFunctionBegin;
18488aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
18588aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
1868c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1872af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
18888aca1feSMatthew G. Knepley   if (!minHeight && dmAux) {
18988aca1feSMatthew G. Knepley     DMLabel spmap;
19088aca1feSMatthew G. Knepley 
19188aca1feSMatthew G. Knepley     /* If dmAux is a surface, then force the projection to take place over a surface */
19288aca1feSMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
19388aca1feSMatthew G. Knepley     if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;}
19488aca1feSMatthew G. Knepley   }
1958c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
1962af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
1978c6c5593SMatthew G. Knepley   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %d not in [0, %d)\n", maxHeight, dim);}
1980c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1998c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
20047923291SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
20147923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2020c898c18SMatthew G. Knepley   if (dmAux) {
2030c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
2040c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
2050c898c18SMatthew G. Knepley   }
206496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
207496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
2088c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
2098c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
2108c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
2118c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
21247923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
21347923291SMatthew G. Knepley     PetscObject  obj;
21447923291SMatthew G. Knepley     PetscClassId id;
21547923291SMatthew G. Knepley 
21647923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
21747923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
21847923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
21947923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
22047923291SMatthew G. Knepley 
22147923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
22247923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
22347923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
22447923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
22547923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
22647923291SMatthew G. Knepley 
22747923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
22847923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
22947923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
23047923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
23147923291SMatthew G. Knepley   }
2320c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
23374fc349aSMatthew G. Knepley     PetscInt   effectiveHeight = auxBd ? minHeight : 0;
23474fc349aSMatthew G. Knepley     PetscFE    fem, subfem;
2350c898c18SMatthew G. Knepley     PetscReal *points;
23674fc349aSMatthew G. Knepley     PetscInt   numPoints, spDim, qdim = 0, d;
2370c898c18SMatthew G. Knepley 
2382af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
2390c898c18SMatthew G. Knepley     numPoints = 0;
2400c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
24174fc349aSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
24274fc349aSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
24374fc349aSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
2440c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2450c898c18SMatthew G. Knepley         if (funcs[f]) {
2460c898c18SMatthew G. Knepley           PetscQuadrature quad;
2470c898c18SMatthew G. Knepley           PetscInt        Nq;
2480c898c18SMatthew G. Knepley 
24974fc349aSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
25074fc349aSMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
2510c898c18SMatthew G. Knepley           numPoints += Nq;
2520c898c18SMatthew G. Knepley         }
2530c898c18SMatthew G. Knepley       }
2540c898c18SMatthew G. Knepley     }
25574fc349aSMatthew G. Knepley     ierr = PetscMalloc1(numPoints*qdim, &points);CHKERRQ(ierr);
2560c898c18SMatthew G. Knepley     numPoints = 0;
2570c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
25874fc349aSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
2590c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2600c898c18SMatthew G. Knepley         if (funcs[f]) {
2610c898c18SMatthew G. Knepley           PetscQuadrature  quad;
2620c898c18SMatthew G. Knepley           const PetscReal *qpoints;
2630c898c18SMatthew G. Knepley           PetscInt         Nq, q;
2640c898c18SMatthew G. Knepley 
26574fc349aSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
2669c3cf19fSMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr);
26774fc349aSMatthew G. Knepley           for (q = 0; q < Nq*qdim; ++q) points[numPoints*qdim+q] = qpoints[q];
2680c898c18SMatthew G. Knepley           numPoints += Nq;
2690c898c18SMatthew G. Knepley         }
2700c898c18SMatthew G. Knepley       }
2710c898c18SMatthew G. Knepley     }
2720c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
2730c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
2740c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
2750c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
27674fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
27774fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
27874fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
2790c898c18SMatthew G. Knepley     }
2800c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
2810c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
28274fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
28374fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
28474fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
2850c898c18SMatthew G. Knepley     }
2860c898c18SMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
2870c898c18SMatthew G. Knepley   }
28847923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
2892af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
29088aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
2918c6c5593SMatthew G. Knepley     PetscScalar *values;
2928c6c5593SMatthew G. Knepley     PetscBool   *fieldActive;
29384f7fa7aSMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
2948c6c5593SMatthew G. Knepley 
29547923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
2968c6c5593SMatthew G. Knepley     if (!h) {
2978c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
2988c6c5593SMatthew G. Knepley 
2998c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3008c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
3018c6c5593SMatthew G. Knepley     }
30247923291SMatthew G. Knepley     if (pEnd <= pStart) continue;
30347923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
30447923291SMatthew G. Knepley     totDim = 0;
30547923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
306bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
30747923291SMatthew G. Knepley         sp[f] = cellsp[f];
30847923291SMatthew G. Knepley       } else {
309bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
31047923291SMatthew G. Knepley         if (!sp[f]) continue;
31147923291SMatthew G. Knepley       }
31247923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
3139c3cf19fSMatthew G. Knepley       totDim += spDim;
31447923291SMatthew G. Knepley     }
31547923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
3168c6c5593SMatthew G. Knepley     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
31747923291SMatthew G. Knepley     if (!totDim) continue;
31847923291SMatthew G. Knepley     /* Loop over points at this height */
319*69291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
320*69291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
3218c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
3228c6c5593SMatthew G. Knepley     if (label) {
3238c6c5593SMatthew G. Knepley       PetscInt i;
32447923291SMatthew G. Knepley 
32547923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
32647923291SMatthew G. Knepley         IS              pointIS;
32747923291SMatthew G. Knepley         const PetscInt *points;
3288c6c5593SMatthew G. Knepley         PetscInt        n;
32947923291SMatthew G. Knepley 
33047923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
33147923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
33247923291SMatthew G. Knepley         ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
33347923291SMatthew G. Knepley         ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
33447923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
33547923291SMatthew G. Knepley           const PetscInt  point = points[p];
33647923291SMatthew G. Knepley 
33747923291SMatthew G. Knepley           if ((point < pStart) || (point >= pEnd)) continue;
338b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
3391c531cf8SMatthew G. Knepley           ierr = DMProjectPoint_Private(dm, dmAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
34047923291SMatthew G. Knepley           if (ierr) {
34147923291SMatthew G. Knepley             PetscErrorCode ierr2;
342*69291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
343*69291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
34447923291SMatthew G. Knepley             CHKERRQ(ierr);
34547923291SMatthew G. Knepley           }
34647923291SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
34747923291SMatthew G. Knepley         }
34847923291SMatthew G. Knepley         ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
34947923291SMatthew G. Knepley         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
35047923291SMatthew G. Knepley       }
3518c6c5593SMatthew G. Knepley     } else {
3528c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
353b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
3541c531cf8SMatthew G. Knepley         ierr = DMProjectPoint_Private(dm, dmAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
3558c6c5593SMatthew G. Knepley         if (ierr) {
3568c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
357*69291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
358*69291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
3598c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
3608c6c5593SMatthew G. Knepley         }
3618c6c5593SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, values, mode);CHKERRQ(ierr);
3628c6c5593SMatthew G. Knepley       }
3638c6c5593SMatthew G. Knepley     }
364*69291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
365*69291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
36647923291SMatthew G. Knepley   }
3678c6c5593SMatthew G. Knepley   /* Cleanup */
3680c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
36974fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
37074fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
3710c898c18SMatthew G. Knepley 
3720c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3730c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3740c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
37574fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
37674fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
37774fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3780c898c18SMatthew G. Knepley     }
3790c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3800c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
38174fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
38274fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
38374fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3840c898c18SMatthew G. Knepley     }
3850c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
3860c898c18SMatthew G. Knepley   }
387496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
3888c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
3898c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
39047923291SMatthew G. Knepley }
3918c6c5593SMatthew G. Knepley 
3928c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
3938c6c5593SMatthew G. Knepley {
3948c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
3958c6c5593SMatthew G. Knepley 
3968c6c5593SMatthew G. Knepley   PetscFunctionBegin;
3971c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
3988c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
3998c6c5593SMatthew G. Knepley }
4008c6c5593SMatthew G. Knepley 
4011c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
4028c6c5593SMatthew G. Knepley {
4038c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
4048c6c5593SMatthew G. Knepley 
4058c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4061c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
40747923291SMatthew G. Knepley   PetscFunctionReturn(0);
40847923291SMatthew G. Knepley }
40947923291SMatthew G. Knepley 
4108c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
41147923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
41247923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
41347923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
414191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
41547923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
41647923291SMatthew G. Knepley {
41747923291SMatthew G. Knepley   PetscErrorCode ierr;
41847923291SMatthew G. Knepley 
41947923291SMatthew G. Knepley   PetscFunctionBegin;
4201c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
4218c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
42247923291SMatthew G. Knepley }
42347923291SMatthew G. Knepley 
4241c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
4258c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
4268c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
4278c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
428191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
4298c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
4308c6c5593SMatthew G. Knepley {
4318c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
43247923291SMatthew G. Knepley 
4338c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4341c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
43547923291SMatthew G. Knepley   PetscFunctionReturn(0);
43647923291SMatthew G. Knepley }
437