1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petsc/private/petscfeimpl.h> 4 5 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS prob, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], 6 PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 7 PetscScalar values[]) 8 { 9 PetscInt coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp; 10 PetscBool isAffine, transform; 11 PetscErrorCode ierr; 12 13 PetscFunctionBeginHot; 14 ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr); 15 ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 16 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 17 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 18 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 19 /* Get values for closure */ 20 isAffine = fegeom->isAffine; 21 for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 22 void * const ctx = ctxs ? ctxs[f] : NULL; 23 24 if (!sp[f]) continue; 25 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 26 if (funcs[f]) { 27 if (isFE[f]) { 28 PetscQuadrature allPoints; 29 PetscInt q, dim, numPoints; 30 const PetscReal *points; 31 PetscScalar *pointEval; 32 PetscReal *x; 33 DM rdm; 34 35 ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr); 36 ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 37 ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 38 ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 39 ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 40 for (q = 0; q < numPoints; q++, tp++) { 41 const PetscReal *v0; 42 43 if (isAffine) { 44 CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x); 45 v0 = x; 46 } else { 47 v0 = &fegeom->v[tp*coordDim]; 48 } 49 if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;} 50 ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 51 } 52 /* Transform point evaluations pointEval[q,c] */ 53 ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr); 54 ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 55 ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 56 ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 57 v += spDim; 58 } else { 59 for (d = 0; d < spDim; ++d, ++v) { 60 ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 61 } 62 } 63 } else { 64 for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;} 65 } 66 } 67 PetscFunctionReturn(0); 68 } 69 70 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 71 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 72 void (**funcs)(PetscInt, PetscInt, PetscInt, 73 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 74 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 75 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 76 PetscScalar values[]) 77 { 78 PetscSection section, sectionAux = NULL; 79 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 80 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 81 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 82 const PetscScalar *constants; 83 PetscReal *x; 84 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 85 PetscFEGeom fegeom; 86 const PetscInt dE = cgeom->dimEmbed; 87 PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 88 PetscBool isAffine, transform; 89 PetscErrorCode ierr; 90 91 PetscFunctionBeginHot; 92 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 93 ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 94 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 95 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 96 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 97 ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 98 ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 99 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 100 ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 101 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 102 ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 103 if (dmAux) { 104 PetscInt subp; 105 106 ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr); 107 ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 108 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 109 ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 110 ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 111 ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 112 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 113 ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 114 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 115 ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 116 } 117 /* Get values for closure */ 118 isAffine = cgeom->isAffine; 119 if (isAffine) { 120 fegeom.v = x; 121 fegeom.xi = cgeom->xi; 122 fegeom.J = cgeom->J; 123 fegeom.invJ = cgeom->invJ; 124 fegeom.detJ = cgeom->detJ; 125 } 126 for (f = 0, v = 0; f < Nf; ++f) { 127 PetscQuadrature allPoints; 128 PetscInt q, dim, numPoints; 129 const PetscReal *points; 130 PetscScalar *pointEval; 131 DM dm; 132 133 if (!sp[f]) continue; 134 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 135 if (!funcs[f]) { 136 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 137 continue; 138 } 139 ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 140 ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 141 ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 142 ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 143 for (q = 0; q < numPoints; ++q, ++tp) { 144 if (isAffine) { 145 CoordinatesRefToReal(dE, dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x); 146 } else { 147 fegeom.v = &cgeom->v[tp*dE]; 148 fegeom.J = &cgeom->J[tp*dE*dE]; 149 fegeom.invJ = &cgeom->invJ[tp*dE*dE]; 150 fegeom.detJ = &cgeom->detJ[tp]; 151 } 152 ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 153 if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 154 if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);} 155 (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f]*q]); 156 } 157 ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 158 ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 159 v += spDim; 160 } 161 ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 162 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 163 PetscFunctionReturn(0); 164 } 165 166 static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 167 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 168 void (**funcs)(PetscInt, PetscInt, PetscInt, 169 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 170 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 171 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 172 PetscScalar values[]) 173 { 174 PetscSection section, sectionAux = NULL; 175 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 176 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 177 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 178 const PetscScalar *constants; 179 PetscReal *x; 180 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 181 PetscFEGeom fegeom, cgeom; 182 const PetscInt dE = fgeom->dimEmbed; 183 PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 184 PetscBool isAffine; 185 PetscErrorCode ierr; 186 187 PetscFunctionBeginHot; 188 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 189 ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 190 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 191 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 192 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 193 ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 194 ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 195 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 196 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 197 ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 198 if (dmAux) { 199 DMLabel spmap; 200 PetscInt subp = p; 201 202 /* If dm is a submesh, do not get subpoint */ 203 ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr); 204 if (!spmap) {ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);} 205 ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 206 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 207 ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 208 ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 209 ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 210 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 211 ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 212 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 213 ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 214 } 215 /* Get values for closure */ 216 isAffine = fgeom->isAffine; 217 if (isAffine) { 218 fegeom.v = x; 219 fegeom.xi = fgeom->xi; 220 fegeom.J = fgeom->J; 221 fegeom.invJ = fgeom->invJ; 222 fegeom.detJ = fgeom->detJ; 223 fegeom.n = fgeom->n; 224 225 cgeom.J = fgeom->suppJ[0]; 226 cgeom.invJ = fgeom->suppInvJ[0]; 227 cgeom.detJ = fgeom->suppDetJ[0]; 228 } 229 for (f = 0, v = 0; f < Nf; ++f) { 230 PetscQuadrature allPoints; 231 PetscInt q, dim, numPoints; 232 const PetscReal *points; 233 PetscScalar *pointEval; 234 DM dm; 235 236 if (!sp[f]) continue; 237 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 238 if (!funcs[f]) { 239 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 240 continue; 241 } 242 ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 243 ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 244 ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 245 ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 246 for (q = 0; q < numPoints; ++q, ++tp) { 247 if (isAffine) { 248 CoordinatesRefToReal(dE, dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x); 249 } else { 250 fegeom.v = &fegeom.v[tp*dE]; 251 fegeom.J = &fgeom->J[tp*dE*dE]; 252 fegeom.invJ = &fgeom->invJ[tp*dE*dE]; 253 fegeom.detJ = &fgeom->detJ[tp]; 254 fegeom.n = &fgeom->n[tp*dE]; 255 256 cgeom.J = &fgeom->suppJ[0][tp*dE*dE]; 257 cgeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE]; 258 cgeom.detJ = &fgeom->suppDetJ[0][tp]; 259 } 260 ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 261 if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 262 (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f]*q]); 263 } 264 ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 265 ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 266 v += spDim; 267 } 268 ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 269 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 270 PetscFunctionReturn(0); 271 } 272 273 static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, PetscFEGeom *fegeom, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], 274 PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 275 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 276 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 277 { 278 PetscFVCellGeom fvgeom; 279 PetscInt dim, dimEmbed; 280 PetscErrorCode ierr; 281 282 PetscFunctionBeginHot; 283 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 284 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 285 if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 286 switch (type) { 287 case DM_BC_ESSENTIAL: 288 case DM_BC_NATURAL: 289 ierr = DMProjectPoint_Func_Private(dm, prob, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break; 290 case DM_BC_ESSENTIAL_FIELD: 291 case DM_BC_NATURAL_FIELD: 292 ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps, 293 basisTab, basisDerTab, basisTabAux, basisDerTabAux, 294 (void (**)(PetscInt, PetscInt, PetscInt, 295 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 296 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 297 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 298 case DM_BC_ESSENTIAL_BD_FIELD: 299 ierr = DMProjectPoint_BdField_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps, 300 basisTab, basisDerTab, basisTabAux, basisDerTabAux, 301 (void (**)(PetscInt, PetscInt, PetscInt, 302 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 303 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 304 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 305 default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 306 } 307 PetscFunctionReturn(0); 308 } 309 310 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 311 { 312 PetscReal *points; 313 PetscInt f, numPoints; 314 PetscErrorCode ierr; 315 316 PetscFunctionBegin; 317 numPoints = 0; 318 for (f = 0; f < Nf; ++f) { 319 if (funcs[f]) { 320 PetscQuadrature fAllPoints; 321 PetscInt fNumPoints; 322 323 ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 324 ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 325 numPoints += fNumPoints; 326 } 327 } 328 ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 329 numPoints = 0; 330 for (f = 0; f < Nf; ++f) { 331 PetscInt spDim; 332 333 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 334 if (funcs[f]) { 335 PetscQuadrature fAllPoints; 336 PetscInt qdim, fNumPoints, q; 337 const PetscReal *fPoints; 338 339 ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 340 ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 341 if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim); 342 for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 343 numPoints += fNumPoints; 344 } 345 } 346 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 347 ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 348 PetscFunctionReturn(0); 349 } 350 351 static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob) 352 { 353 DMLabel depthLabel; 354 PetscInt dim, cdepth, ls = -1, i; 355 PetscErrorCode ierr; 356 357 PetscFunctionBegin; 358 if (lStart) *lStart = -1; 359 if (!label) PetscFunctionReturn(0); 360 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 361 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 362 cdepth = dim - height; 363 for (i = 0; i < numIds; ++i) { 364 IS pointIS; 365 const PetscInt *points; 366 PetscInt pdepth; 367 368 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 369 if (!pointIS) continue; /* No points with that id on this process */ 370 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 371 ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr); 372 if (pdepth == cdepth) { 373 ls = points[0]; 374 if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);} 375 } 376 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 377 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 378 if (ls >= 0) break; 379 } 380 if (lStart) *lStart = ls; 381 PetscFunctionReturn(0); 382 } 383 384 /* 385 This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 386 387 There are several different scenarios: 388 389 1) Volumetric mesh with volumetric auxiliary data 390 391 Here minHeight=0 since we loop over cells. 392 393 2) Boundary mesh with boundary auxiliary data 394 395 Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation. 396 397 3) Volumetric mesh with boundary auxiliary data 398 399 Here minHeight=1 and auxbd=PETSC_TRUE since we loop over faces and use data only supported on those faces. This is common when imposing Dirichlet boundary conditions. 400 401 The maxHeight is used to support enforcement of constraints in DMForest. 402 403 If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 404 405 If we are using an input field (DM_BC_ESSENTIAL_FIELD or DM_BC_NATURAL_FIELD), we need to evaluate it at all the quadrature points of the dual basis functionals. 406 - We use effectiveHeight to mean the height above our incoming DS. For example, if the DS is for a submesh then the effective height is zero, whereas if the DS 407 is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract 408 dual spaces for the boundary from our input spaces. 409 - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 410 411 We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 412 413 If we have a label, we iterate over those points. This will probably break the maxHeight functionality since we do not check the height of those points. 414 */ 415 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 416 PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 417 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 418 InsertMode mode, Vec localX) 419 { 420 DM dmAux = NULL, tdm; 421 PetscDS prob = NULL, probAux = NULL; 422 Vec localA = NULL, tv; 423 PetscSection section; 424 PetscDualSpace *sp, *cellsp; 425 PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 426 PetscInt *Nc; 427 PetscInt dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f; 428 PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform; 429 DMField coordField; 430 DMLabel depthLabel; 431 PetscQuadrature allPoints = NULL; 432 PetscErrorCode ierr; 433 434 PetscFunctionBegin; 435 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 436 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 437 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 438 ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr); 439 ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 440 ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 441 ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 442 /* Auxiliary information can only be used with interpolation of field functions */ 443 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 444 if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector"); 445 if (!minHeight && dmAux) { 446 DMLabel spmap; 447 448 /* If dmAux is a surface, then force the projection to take place over a surface */ 449 ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr); 450 if (spmap) { 451 ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr); 452 auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE; 453 } 454 } 455 } 456 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 457 ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr); 458 ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr); 459 maxHeight = PetscMax(maxHeight, minHeight); 460 if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 461 ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr); 462 if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);} 463 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 464 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 465 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 466 if (dmAux) { 467 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 468 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 469 } 470 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 471 ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 472 if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 473 else {cellsp = sp;} 474 if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 475 /* Get cell dual spaces */ 476 for (f = 0; f < Nf; ++f) { 477 PetscObject obj; 478 PetscClassId id; 479 480 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 481 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 482 if (id == PETSCFE_CLASSID) { 483 PetscFE fe = (PetscFE) obj; 484 485 hasFE = PETSC_TRUE; 486 isFE[f] = PETSC_TRUE; 487 ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 488 } else if (id == PETSCFV_CLASSID) { 489 PetscFV fv = (PetscFV) obj; 490 491 hasFV = PETSC_TRUE; 492 isFE[f] = PETSC_FALSE; 493 ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 494 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 495 } 496 ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 497 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 498 PetscInt effectiveHeight = auxBd ? minHeight : 0; 499 PetscFE fem, subfem; 500 const PetscReal *points; 501 PetscInt numPoints; 502 503 if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 504 for (f = 0; f < Nf; ++f) { 505 if (!effectiveHeight) {sp[f] = cellsp[f];} 506 else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 507 } 508 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr); 509 ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 510 ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 511 for (f = 0; f < Nf; ++f) { 512 if (!isFE[f]) continue; 513 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 514 if (!effectiveHeight) {subfem = fem;} 515 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 516 ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 517 } 518 for (f = 0; f < NfAux; ++f) { 519 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 520 if (!effectiveHeight || auxBd) {subfem = fem;} 521 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 522 ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 523 } 524 } 525 /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 526 for (h = minHeight; h <= maxHeight; h++) { 527 PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 528 PetscDS probEff = prob; 529 PetscScalar *values; 530 PetscBool *fieldActive; 531 PetscInt maxDegree; 532 PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 533 IS heightIS; 534 535 ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 536 ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 537 ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 538 if (!h) { 539 PetscInt cEndInterior; 540 541 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 542 pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 543 } 544 if (pEnd <= pStart) { 545 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 546 continue; 547 } 548 /* Compute totDim, the number of dofs in the closure of a point at this height */ 549 totDim = 0; 550 for (f = 0; f < Nf; ++f) { 551 if (!effectiveHeight) { 552 sp[f] = cellsp[f]; 553 } else { 554 ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 555 if (!sp[f]) continue; 556 } 557 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 558 totDim += spDim; 559 } 560 ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr); 561 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim); 562 if (!totDim) { 563 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 564 continue; 565 } 566 if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);} 567 /* Loop over points at this height */ 568 ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 569 ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 570 for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 571 if (label) { 572 PetscInt i; 573 574 for (i = 0; i < numIds; ++i) { 575 IS pointIS, isectIS; 576 const PetscInt *points; 577 PetscInt n; 578 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 579 PetscQuadrature quad = NULL; 580 581 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 582 if (!pointIS) continue; /* No points with that id on this process */ 583 ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 584 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 585 if (!isectIS) continue; 586 ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 587 ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 588 ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 589 if (maxDegree <= 1) { 590 ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 591 } 592 if (!quad) { 593 if (!h && allPoints) { 594 quad = allPoints; 595 allPoints = NULL; 596 } else { 597 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 598 } 599 } 600 ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 601 for (p = 0; p < n; ++p) { 602 const PetscInt point = points[p]; 603 604 ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 605 ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 606 ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 607 if (ierr) { 608 PetscErrorCode ierr2; 609 ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 610 ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 611 CHKERRQ(ierr); 612 } 613 if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 614 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr); 615 } 616 ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 617 ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 618 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 619 ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 620 ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 621 } 622 } else { 623 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 624 PetscQuadrature quad = NULL; 625 IS pointIS; 626 627 ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 628 ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 629 if (maxDegree <= 1) { 630 ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 631 } 632 if (!quad) { 633 if (!h && allPoints) { 634 quad = allPoints; 635 allPoints = NULL; 636 } else { 637 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 638 } 639 } 640 ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 641 for (p = pStart; p < pEnd; ++p) { 642 ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 643 ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 644 ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 645 if (ierr) { 646 PetscErrorCode ierr2; 647 ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 648 ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 649 CHKERRQ(ierr); 650 } 651 if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 652 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr); 653 } 654 ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 655 ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 656 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 657 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 658 } 659 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 660 ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 661 ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 662 } 663 /* Cleanup */ 664 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 665 PetscInt effectiveHeight = auxBd ? minHeight : 0; 666 PetscFE fem, subfem; 667 668 for (f = 0; f < Nf; ++f) { 669 if (!isFE[f]) continue; 670 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 671 if (!effectiveHeight) {subfem = fem;} 672 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 673 ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 674 } 675 for (f = 0; f < NfAux; ++f) { 676 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 677 if (!effectiveHeight || auxBd) {subfem = fem;} 678 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 679 ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 680 } 681 ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 682 } 683 ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 684 ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 685 if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 686 PetscFunctionReturn(0); 687 } 688 689 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 690 { 691 PetscErrorCode ierr; 692 693 PetscFunctionBegin; 694 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 695 PetscFunctionReturn(0); 696 } 697 698 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) 699 { 700 PetscErrorCode ierr; 701 702 PetscFunctionBegin; 703 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 704 PetscFunctionReturn(0); 705 } 706 707 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 708 void (**funcs)(PetscInt, PetscInt, PetscInt, 709 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 710 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 711 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 712 InsertMode mode, Vec localX) 713 { 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720 721 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 722 void (**funcs)(PetscInt, PetscInt, PetscInt, 723 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 724 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 725 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 726 InsertMode mode, Vec localX) 727 { 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 732 PetscFunctionReturn(0); 733 } 734