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