xref: /petsc/src/dm/impls/plex/plexproject.c (revision 496733e227f8aa6258a82f10fb5673d2cffcc35b)
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 
547923291SMatthew G. Knepley #undef __FUNCT__
68c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectPoint_Func_Private"
7*496733e2SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscReal time, PetscFECellGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
88c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
98c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
1047923291SMatthew G. Knepley {
118c6c5593SMatthew G. Knepley   PetscDS        prob;
12*496733e2SMatthew G. Knepley   PetscInt       Nf, *Nc, f, spDim, d, c, v;
138c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
148c6c5593SMatthew G. Knepley 
158c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
168c6c5593SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
178c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
18*496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
198c6c5593SMatthew G. Knepley   /* Get values for closure */
208c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
218c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
228c6c5593SMatthew G. Knepley 
238c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
248c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
258c6c5593SMatthew G. Knepley     for (d = 0; d < spDim; ++d) {
268c6c5593SMatthew G. Knepley       if (funcs[f]) {
278c6c5593SMatthew G. Knepley         if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, fegeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
288c6c5593SMatthew G. Knepley         else         {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
298c6c5593SMatthew G. Knepley       } else {
308c6c5593SMatthew G. Knepley         for (c = 0; c < Nc[f]; ++c) values[v+c] = 0.0;
318c6c5593SMatthew G. Knepley       }
328c6c5593SMatthew G. Knepley       v += Nc[f];
338c6c5593SMatthew G. Knepley     }
348c6c5593SMatthew G. Knepley   }
358c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
368c6c5593SMatthew G. Knepley }
378c6c5593SMatthew G. Knepley 
388c6c5593SMatthew G. Knepley #undef __FUNCT__
398c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectPoint_Field_Private"
40*496733e2SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, DM dmAux, PetscReal time, Vec localU, Vec localA, PetscFECellGeom *fegeom, PetscDualSpace sp[], PetscInt p,
410c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
428c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
438c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
458c6c5593SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscScalar[]), void **ctxs,
468c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
478c6c5593SMatthew G. Knepley {
488c6c5593SMatthew G. Knepley   PetscDS        prob, probAux = NULL;
498c6c5593SMatthew G. Knepley   PetscSection   section, sectionAux = NULL;
50*496733e2SMatthew G. Knepley   PetscScalar   *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL;
518c6c5593SMatthew G. Knepley   PetscScalar   *coefficients   = NULL, *coefficientsAux   = NULL;
528c6c5593SMatthew G. Knepley   PetscScalar   *coefficients_t = NULL, *coefficientsAux_t = NULL;
538c6c5593SMatthew G. Knepley   PetscReal     *x;
54*496733e2SMatthew G. Knepley   PetscInt      *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
550c898c18SMatthew G. Knepley   PetscInt       dim, Nf, NfAux = 0, f, spDim, d, c, v, tp = 0;
568c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
578c6c5593SMatthew G. Knepley 
588c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
598c6c5593SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
608c6c5593SMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
618c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
62*496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
63*496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
648c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
658c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
6611d189daSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL /*&u_t*/, &u_x);CHKERRQ(ierr);
67*496733e2SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
688c6c5593SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
69*496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
708c6c5593SMatthew G. Knepley   if (dmAux) {
718c6c5593SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
728c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
73*496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
74*496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
758c6c5593SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
768c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
778c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
7811d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
79*496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
80*496733e2SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);
818c6c5593SMatthew G. Knepley   }
828c6c5593SMatthew G. Knepley   /* Get values for closure */
838c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
848c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
858c6c5593SMatthew G. Knepley 
868c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
878c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
888c6c5593SMatthew G. Knepley     for (d = 0; d < spDim; ++d) {
898c6c5593SMatthew G. Knepley       if (funcs[f]) {
908c6c5593SMatthew G. Knepley         PetscQuadrature  quad;
910c898c18SMatthew G. Knepley         const PetscReal *points;
928c6c5593SMatthew G. Knepley         PetscInt         numPoints, q;
938c6c5593SMatthew G. Knepley 
948c6c5593SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
950c898c18SMatthew G. Knepley         ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, NULL);CHKERRQ(ierr);
960c898c18SMatthew G. Knepley         for (q = 0; q < numPoints; ++q, ++tp) {
978c6c5593SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, fegeom->v0, fegeom->J, &points[q*dim], x);
98*496733e2SMatthew G. Knepley           EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t);
99*496733e2SMatthew G. Knepley           if (probAux) {EvaluateFieldJets(dim, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
1008c6c5593SMatthew 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, &values[v]);
1018c6c5593SMatthew G. Knepley         }
1028c6c5593SMatthew G. Knepley       } else {
1038c6c5593SMatthew G. Knepley         for (c = 0; c < Nc[f]; ++c) values[v+c] = 0.0;
1048c6c5593SMatthew G. Knepley       }
1058c6c5593SMatthew G. Knepley       v += Nc[f];
1068c6c5593SMatthew G. Knepley     }
1078c6c5593SMatthew G. Knepley   }
1088c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1098c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1108c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1118c6c5593SMatthew G. Knepley }
1128c6c5593SMatthew G. Knepley 
1138c6c5593SMatthew G. Knepley #undef __FUNCT__
1148c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectPoint_Private"
1150c898c18SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Private(DM dm, DM dmAux, PetscInt h, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
116*496733e2SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1178c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1188c6c5593SMatthew G. Knepley {
1198c6c5593SMatthew G. Knepley   PetscFECellGeom fegeom;
1208c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1218c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1228c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1238c6c5593SMatthew G. Knepley 
1248c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1258c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1268c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1278c6c5593SMatthew G. Knepley   if (hasFE) {
12843048b12SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr);
1298c6c5593SMatthew G. Knepley     fegeom.dim      = dim - h;
1308c6c5593SMatthew G. Knepley     fegeom.dimEmbed = dimEmbed;
1318c6c5593SMatthew G. Knepley   }
1328c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1338c6c5593SMatthew G. Knepley   switch (type) {
1348c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1358c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
136*496733e2SMatthew 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;
1378c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1388c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
139*496733e2SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, dmAux, time, localU, localA, &fegeom, sp, p,
1400c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1410c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1420c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1430c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1440c898c18SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1458c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1468c6c5593SMatthew G. Knepley   }
1478c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1488c6c5593SMatthew G. Knepley }
1498c6c5593SMatthew G. Knepley 
1508c6c5593SMatthew G. Knepley #undef __FUNCT__
1518c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectLocal_Generic_Plex"
1528c6c5593SMatthew G. Knepley PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
1538c6c5593SMatthew G. Knepley                                            DMLabel label, PetscInt numIds, const PetscInt ids[],
1548c6c5593SMatthew G. Knepley                                            DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
1558c6c5593SMatthew G. Knepley                                            InsertMode mode, Vec localX)
1568c6c5593SMatthew G. Knepley {
1578c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
1580c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
1598c6c5593SMatthew G. Knepley   Vec             localA = NULL;
16047923291SMatthew G. Knepley   PetscSection    section;
1618c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
1620c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
1638c6c5593SMatthew G. Knepley   PetscInt       *Nc;
1640c898c18SMatthew G. Knepley   PetscInt        dim, dimEmbed, maxHeight, h, Nf, NfAux = 0, f;
16547923291SMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE;
16647923291SMatthew G. Knepley   PetscErrorCode  ierr;
16747923291SMatthew G. Knepley 
16847923291SMatthew G. Knepley   PetscFunctionBegin;
1698c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1708c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
1718c6c5593SMatthew 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);}
1720c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1738c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
17447923291SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
17547923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1768c6c5593SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1778c6c5593SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
1780c898c18SMatthew G. Knepley   if (dmAux) {
1790c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1800c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
1810c898c18SMatthew G. Knepley   }
182*496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
183*496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
1848c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
1858c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
1868c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
1878c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
18847923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
18947923291SMatthew G. Knepley     PetscObject  obj;
19047923291SMatthew G. Knepley     PetscClassId id;
19147923291SMatthew G. Knepley 
19247923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
19347923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
19447923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
19547923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
19647923291SMatthew G. Knepley 
19747923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
19847923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
19947923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
20047923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
20147923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
20247923291SMatthew G. Knepley 
20347923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
20447923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
20547923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
20647923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
20747923291SMatthew G. Knepley   }
2080c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
2090c898c18SMatthew G. Knepley     PetscFE    fem;
2100c898c18SMatthew G. Knepley     PetscReal *points;
2110c898c18SMatthew G. Knepley     PetscInt   numPoints, spDim, d;
2120c898c18SMatthew G. Knepley 
2130c898c18SMatthew G. Knepley     if (maxHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field proejction not supported for face interpolation");
2140c898c18SMatthew G. Knepley     numPoints = 0;
2150c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
2160c898c18SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr);
2170c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2180c898c18SMatthew G. Knepley         if (funcs[f]) {
2190c898c18SMatthew G. Knepley           PetscQuadrature quad;
2200c898c18SMatthew G. Knepley           PetscInt        Nq;
2210c898c18SMatthew G. Knepley 
2220c898c18SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr);
2230c898c18SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
2240c898c18SMatthew G. Knepley           numPoints += Nq;
2250c898c18SMatthew G. Knepley         }
2260c898c18SMatthew G. Knepley       }
2270c898c18SMatthew G. Knepley     }
2280c898c18SMatthew G. Knepley     ierr = PetscMalloc1(numPoints*dim, &points);CHKERRQ(ierr);
2290c898c18SMatthew G. Knepley     numPoints = 0;
2300c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
2310c898c18SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr);
2320c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2330c898c18SMatthew G. Knepley         if (funcs[f]) {
2340c898c18SMatthew G. Knepley           PetscQuadrature  quad;
2350c898c18SMatthew G. Knepley           const PetscReal *qpoints;
2360c898c18SMatthew G. Knepley           PetscInt         Nq, q;
2370c898c18SMatthew G. Knepley 
2380c898c18SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr);
2390c898c18SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr);
2400c898c18SMatthew G. Knepley           for (q = 0; q < Nq*dim; ++q) points[numPoints*dim+q] = qpoints[q];
2410c898c18SMatthew G. Knepley           numPoints += Nq;
2420c898c18SMatthew G. Knepley         }
2430c898c18SMatthew G. Knepley       }
2440c898c18SMatthew G. Knepley     }
2450c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
2460c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
2470c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
2480c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
2490c898c18SMatthew G. Knepley       ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
2500c898c18SMatthew G. Knepley     }
2510c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
2520c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
2530c898c18SMatthew G. Knepley       ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
2540c898c18SMatthew G. Knepley     }
2550c898c18SMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
2560c898c18SMatthew G. Knepley   }
25747923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
25847923291SMatthew G. Knepley   for (h = 0; h <= maxHeight; h++) {
2598c6c5593SMatthew G. Knepley     PetscScalar *values;
2608c6c5593SMatthew G. Knepley     PetscBool   *fieldActive;
2618c6c5593SMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, d, numValues, v;
2628c6c5593SMatthew G. Knepley 
26347923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
2648c6c5593SMatthew G. Knepley     if (!h) {
2658c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
2668c6c5593SMatthew G. Knepley 
2678c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2688c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
2698c6c5593SMatthew G. Knepley     }
27047923291SMatthew G. Knepley     if (pEnd <= pStart) continue;
27147923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
27247923291SMatthew G. Knepley     totDim = 0;
27347923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
27447923291SMatthew G. Knepley       if (!h) {
27547923291SMatthew G. Knepley         sp[f] = cellsp[f];
27647923291SMatthew G. Knepley       } else {
27747923291SMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
27847923291SMatthew G. Knepley         if (!sp[f]) continue;
27947923291SMatthew G. Knepley       }
28047923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
2818c6c5593SMatthew G. Knepley       totDim += spDim*Nc[f];
28247923291SMatthew G. Knepley     }
28347923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
2848c6c5593SMatthew 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);
28547923291SMatthew G. Knepley     if (!totDim) continue;
28647923291SMatthew G. Knepley     /* Loop over points at this height */
28747923291SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
2888c6c5593SMatthew G. Knepley     ierr = DMGetWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
2898c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
2908c6c5593SMatthew G. Knepley     if (label) {
2918c6c5593SMatthew G. Knepley       PetscInt i;
29247923291SMatthew G. Knepley 
29347923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
29447923291SMatthew G. Knepley         IS              pointIS;
29547923291SMatthew G. Knepley         const PetscInt *points;
2968c6c5593SMatthew G. Knepley         PetscInt        n;
29747923291SMatthew G. Knepley 
29847923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
29947923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
30047923291SMatthew G. Knepley         ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
30147923291SMatthew G. Knepley         ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
30247923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
30347923291SMatthew G. Knepley           const PetscInt  point = points[p];
30447923291SMatthew G. Knepley 
30547923291SMatthew G. Knepley           if ((point < pStart) || (point >= pEnd)) continue;
306*496733e2SMatthew G. Knepley           ierr = DMProjectPoint_Private(dm, dmAux, h, time, localU, localA, hasFE, hasFV, isFE, sp, point, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
30747923291SMatthew G. Knepley           if (ierr) {
30847923291SMatthew G. Knepley             PetscErrorCode ierr2;
30947923291SMatthew G. Knepley             ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
3108c6c5593SMatthew G. Knepley             ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2);
31147923291SMatthew G. Knepley             CHKERRQ(ierr);
31247923291SMatthew G. Knepley           }
31347923291SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
31447923291SMatthew G. Knepley         }
31547923291SMatthew G. Knepley         ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
31647923291SMatthew G. Knepley         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
31747923291SMatthew G. Knepley       }
3188c6c5593SMatthew G. Knepley     } else {
3198c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
320*496733e2SMatthew G. Knepley         ierr = DMProjectPoint_Private(dm, dmAux, h, time, localU, localA, hasFE, hasFV, isFE, sp, p, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
3218c6c5593SMatthew G. Knepley         if (ierr) {
3228c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
3238c6c5593SMatthew G. Knepley           ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
3248c6c5593SMatthew G. Knepley           ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2);
3258c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
3268c6c5593SMatthew G. Knepley         }
3278c6c5593SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, values, mode);CHKERRQ(ierr);
3288c6c5593SMatthew G. Knepley       }
3298c6c5593SMatthew G. Knepley     }
33047923291SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3318c6c5593SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
33247923291SMatthew G. Knepley   }
3338c6c5593SMatthew G. Knepley   /* Cleanup */
3340c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
3350c898c18SMatthew G. Knepley     PetscFE fem;
3360c898c18SMatthew G. Knepley 
3370c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3380c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3390c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
3400c898c18SMatthew G. Knepley       ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3410c898c18SMatthew G. Knepley     }
3420c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3430c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
3440c898c18SMatthew G. Knepley       ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3450c898c18SMatthew G. Knepley     }
3460c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
3470c898c18SMatthew G. Knepley   }
348*496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
3498c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
3508c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
35147923291SMatthew G. Knepley }
3528c6c5593SMatthew G. Knepley 
3538c6c5593SMatthew G. Knepley #undef __FUNCT__
3548c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectFunctionLocal_Plex"
3558c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
3568c6c5593SMatthew G. Knepley {
3578c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
3588c6c5593SMatthew G. Knepley 
3598c6c5593SMatthew G. Knepley   PetscFunctionBegin;
3608c6c5593SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
3618c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
3628c6c5593SMatthew G. Knepley }
3638c6c5593SMatthew G. Knepley 
3648c6c5593SMatthew G. Knepley #undef __FUNCT__
3658c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectFunctionLabelLocal_Plex"
3668c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
3678c6c5593SMatthew G. Knepley {
3688c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
3698c6c5593SMatthew G. Knepley 
3708c6c5593SMatthew G. Knepley   PetscFunctionBegin;
3718c6c5593SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
37247923291SMatthew G. Knepley   PetscFunctionReturn(0);
37347923291SMatthew G. Knepley }
37447923291SMatthew G. Knepley 
37547923291SMatthew G. Knepley #undef __FUNCT__
37647923291SMatthew G. Knepley #define __FUNCT__ "DMProjectFieldLocal_Plex"
3778c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
37847923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
37947923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
38047923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
38147923291SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscScalar[]),
38247923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
38347923291SMatthew G. Knepley {
38447923291SMatthew G. Knepley   PetscErrorCode ierr;
38547923291SMatthew G. Knepley 
38647923291SMatthew G. Knepley   PetscFunctionBegin;
3878c6c5593SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
3888c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
38947923291SMatthew G. Knepley }
39047923291SMatthew G. Knepley 
3918c6c5593SMatthew G. Knepley #undef __FUNCT__
3928c6c5593SMatthew G. Knepley #define __FUNCT__ "DMProjectFieldLabelLocal_Plex"
3938c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], Vec localU,
3948c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
3958c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3968c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3978c6c5593SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscScalar[]),
3988c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
3998c6c5593SMatthew G. Knepley {
4008c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
40147923291SMatthew G. Knepley 
4028c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4038c6c5593SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
40447923291SMatthew G. Knepley   PetscFunctionReturn(0);
40547923291SMatthew G. Knepley }
406