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