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