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, PetscReal time, PetscFECellGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], 6 PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 7 PetscScalar values[]) 8 { 9 PetscDS prob; 10 PetscInt Nf, *Nc, f, spDim, d, c, v; 11 PetscErrorCode ierr; 12 13 PetscFunctionBeginHot; 14 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 15 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 16 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 17 /* Get values for closure */ 18 for (f = 0, v = 0; f < Nf; ++f) { 19 void * const ctx = ctxs ? ctxs[f] : NULL; 20 21 if (!sp[f]) continue; 22 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 23 for (d = 0; d < spDim; ++d) { 24 if (funcs[f]) { 25 if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, fegeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);} 26 else {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);} 27 } else { 28 for (c = 0; c < Nc[f]; ++c) values[v+c] = 0.0; 29 } 30 v += Nc[f]; 31 } 32 } 33 PetscFunctionReturn(0); 34 } 35 36 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, DM dmAux, PetscReal time, Vec localU, Vec localA, PetscFECellGeom *fegeom, PetscDualSpace sp[], PetscInt p, 37 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 38 void (**funcs)(PetscInt, PetscInt, PetscInt, 39 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 40 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 41 PetscReal, const PetscReal[], PetscScalar[]), void **ctxs, 42 PetscScalar values[]) 43 { 44 PetscDS prob, probAux = NULL; 45 PetscSection section, sectionAux = NULL; 46 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL; 47 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 48 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 49 PetscReal *x; 50 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 51 PetscInt dim, Nf, NfAux = 0, f, spDim, d, c, v, tp = 0; 52 PetscErrorCode ierr; 53 54 PetscFunctionBeginHot; 55 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 56 ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 57 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 58 ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 59 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 60 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 61 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 62 ierr = PetscDSGetEvaluationArrays(prob, &u, NULL /*&u_t*/, &u_x);CHKERRQ(ierr); 63 ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 64 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 65 ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 66 if (dmAux) { 67 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 68 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 69 ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 70 ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 71 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 72 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 73 ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 74 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 75 ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 76 ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr); 77 } 78 /* Get values for closure */ 79 for (f = 0, v = 0; f < Nf; ++f) { 80 if (!sp[f]) continue; 81 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 82 for (d = 0; d < spDim; ++d) { 83 if (funcs[f]) { 84 PetscQuadrature quad; 85 const PetscReal *points; 86 PetscInt numPoints, q; 87 88 ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr); 89 ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, NULL);CHKERRQ(ierr); 90 for (q = 0; q < numPoints; ++q, ++tp) { 91 CoordinatesRefToReal(dim, dim, fegeom->v0, fegeom->J, &points[q*dim], x); 92 EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t); 93 if (probAux) {EvaluateFieldJets(dim, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);} 94 (*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]); 95 } 96 } else { 97 for (c = 0; c < Nc[f]; ++c) values[v+c] = 0.0; 98 } 99 v += Nc[f]; 100 } 101 } 102 ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 103 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 104 PetscFunctionReturn(0); 105 } 106 107 static PetscErrorCode DMProjectPoint_Private(DM dm, DM dmAux, PetscInt h, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], 108 PetscDualSpace sp[], PetscInt p, PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 109 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 110 { 111 PetscFECellGeom fegeom; 112 PetscFVCellGeom fvgeom; 113 PetscInt dim, dimEmbed; 114 PetscErrorCode ierr; 115 116 PetscFunctionBeginHot; 117 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 118 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 119 if (hasFE) { 120 ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr); 121 fegeom.dim = dim - h; 122 fegeom.dimEmbed = dimEmbed; 123 } 124 if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 125 switch (type) { 126 case DM_BC_ESSENTIAL: 127 case DM_BC_NATURAL: 128 ierr = DMProjectPoint_Func_Private(dm, time, &fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break; 129 case DM_BC_ESSENTIAL_FIELD: 130 case DM_BC_NATURAL_FIELD: 131 ierr = DMProjectPoint_Field_Private(dm, dmAux, time, localU, localA, &fegeom, sp, p, 132 basisTab, basisDerTab, basisTabAux, basisDerTabAux, 133 (void (**)(PetscInt, PetscInt, PetscInt, 134 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 135 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 136 PetscReal, const PetscReal[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 137 default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 138 } 139 PetscFunctionReturn(0); 140 } 141 142 PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 143 DMLabel label, PetscInt numIds, const PetscInt ids[], 144 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 145 InsertMode mode, Vec localX) 146 { 147 DM dmAux = NULL; 148 PetscDS prob, probAux = NULL; 149 Vec localA = NULL; 150 PetscSection section; 151 PetscDualSpace *sp, *cellsp; 152 PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 153 PetscInt *Nc; 154 PetscInt dim, dimEmbed, minHeight, maxHeight, h, Nf, NfAux = 0, f; 155 PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE; 156 PetscErrorCode ierr; 157 158 PetscFunctionBegin; 159 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 160 ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr); 161 ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr); 162 maxHeight = PetscMax(maxHeight, minHeight); 163 if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %d not in [0, %d)\n", maxHeight, dim);} 164 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 165 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 166 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 167 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 168 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 169 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 170 if (dmAux) { 171 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 172 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 173 } 174 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 175 ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 176 if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 177 else {cellsp = sp;} 178 if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 179 /* Get cell dual spaces */ 180 for (f = 0; f < Nf; ++f) { 181 PetscObject obj; 182 PetscClassId id; 183 184 ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 185 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 186 if (id == PETSCFE_CLASSID) { 187 PetscFE fe = (PetscFE) obj; 188 189 hasFE = PETSC_TRUE; 190 isFE[f] = PETSC_TRUE; 191 ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 192 } else if (id == PETSCFV_CLASSID) { 193 PetscFV fv = (PetscFV) obj; 194 195 hasFV = PETSC_TRUE; 196 isFE[f] = PETSC_FALSE; 197 ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 198 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 199 } 200 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 201 PetscFE fem; 202 PetscReal *points; 203 PetscInt numPoints, spDim, d; 204 205 if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 206 numPoints = 0; 207 for (f = 0; f < Nf; ++f) { 208 ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr); 209 for (d = 0; d < spDim; ++d) { 210 if (funcs[f]) { 211 PetscQuadrature quad; 212 PetscInt Nq; 213 214 ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr); 215 ierr = PetscQuadratureGetData(quad, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 216 numPoints += Nq; 217 } 218 } 219 } 220 ierr = PetscMalloc1(numPoints*dim, &points);CHKERRQ(ierr); 221 numPoints = 0; 222 for (f = 0; f < Nf; ++f) { 223 ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr); 224 for (d = 0; d < spDim; ++d) { 225 if (funcs[f]) { 226 PetscQuadrature quad; 227 const PetscReal *qpoints; 228 PetscInt Nq, q; 229 230 ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr); 231 ierr = PetscQuadratureGetData(quad, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr); 232 for (q = 0; q < Nq*dim; ++q) points[numPoints*dim+q] = qpoints[q]; 233 numPoints += Nq; 234 } 235 } 236 } 237 ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 238 for (f = 0; f < Nf; ++f) { 239 if (!isFE[f]) continue; 240 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 241 ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 242 } 243 for (f = 0; f < NfAux; ++f) { 244 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 245 ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 246 } 247 ierr = PetscFree(points);CHKERRQ(ierr); 248 } 249 /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 250 for (h = minHeight; h <= maxHeight; h++) { 251 PetscScalar *values; 252 PetscBool *fieldActive; 253 PetscInt pStart, pEnd, p, spDim, totDim, numValues; 254 255 ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 256 if (!h) { 257 PetscInt cEndInterior; 258 259 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 260 pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 261 } 262 if (pEnd <= pStart) continue; 263 /* Compute totDim, the number of dofs in the closure of a point at this height */ 264 totDim = 0; 265 for (f = 0; f < Nf; ++f) { 266 PetscInt effectiveHeight = h - minHeight; 267 if (!effectiveHeight) { 268 sp[f] = cellsp[f]; 269 } else { 270 ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 271 if (!sp[f]) continue; 272 } 273 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 274 totDim += spDim*Nc[f]; 275 } 276 ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr); 277 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim); 278 if (!totDim) continue; 279 /* Loop over points at this height */ 280 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 281 ierr = DMGetWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 282 for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 283 if (label) { 284 PetscInt i; 285 286 for (i = 0; i < numIds; ++i) { 287 IS pointIS; 288 const PetscInt *points; 289 PetscInt n; 290 291 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 292 if (!pointIS) continue; /* No points with that id on this process */ 293 ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 294 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 295 for (p = 0; p < n; ++p) { 296 const PetscInt point = points[p]; 297 298 if ((point < pStart) || (point >= pEnd)) continue; 299 ierr = DMProjectPoint_Private(dm, dmAux, h, time, localU, localA, hasFE, hasFV, isFE, sp, point, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 300 if (ierr) { 301 PetscErrorCode ierr2; 302 ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2); 303 ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2); 304 CHKERRQ(ierr); 305 } 306 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 307 } 308 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 309 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 310 } 311 } else { 312 for (p = pStart; p < pEnd; ++p) { 313 ierr = DMProjectPoint_Private(dm, dmAux, h, time, localU, localA, hasFE, hasFV, isFE, sp, p, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 314 if (ierr) { 315 PetscErrorCode ierr2; 316 ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2); 317 ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2); 318 CHKERRQ(ierr); 319 } 320 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, values, mode);CHKERRQ(ierr); 321 } 322 } 323 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 324 ierr = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 325 } 326 /* Cleanup */ 327 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 328 PetscFE fem; 329 330 for (f = 0; f < Nf; ++f) { 331 if (!isFE[f]) continue; 332 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 333 ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 334 } 335 for (f = 0; f < NfAux; ++f) { 336 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 337 ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 338 } 339 ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 340 } 341 ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 342 if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 343 PetscFunctionReturn(0); 344 } 345 346 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 347 { 348 PetscErrorCode ierr; 349 350 PetscFunctionBegin; 351 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 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) 356 { 357 PetscErrorCode ierr; 358 359 PetscFunctionBegin; 360 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 365 void (**funcs)(PetscInt, PetscInt, PetscInt, 366 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 367 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 368 PetscReal, const PetscReal[], PetscScalar[]), 369 InsertMode mode, Vec localX) 370 { 371 PetscErrorCode ierr; 372 373 PetscFunctionBegin; 374 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], Vec localU, 379 void (**funcs)(PetscInt, PetscInt, PetscInt, 380 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 381 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 382 PetscReal, const PetscReal[], PetscScalar[]), 383 InsertMode mode, Vec localX) 384 { 385 PetscErrorCode ierr; 386 387 PetscFunctionBegin; 388 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391