1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petsc/private/petscfeimpl.h> 4 5 /*@ 6 DMPlexGetActivePoint - Get the point on which projection is currently working 7 8 Not Collective 9 10 Input Parameter: 11 . dm - the `DM` 12 13 Output Parameter: 14 . point - The mesh point involved in the current projection 15 16 Level: developer 17 18 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()` 19 @*/ 20 PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point) 21 { 22 PetscFunctionBeginHot; 23 *point = ((DM_Plex *)dm->data)->activePoint; 24 PetscFunctionReturn(PETSC_SUCCESS); 25 } 26 27 /*@ 28 DMPlexSetActivePoint - Set the point on which projection is currently working 29 30 Not Collective 31 32 Input Parameters: 33 + dm - the `DM` 34 - point - The mesh point involved in the current projection 35 36 Level: developer 37 38 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetActivePoint()` 39 @*/ 40 PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point) 41 { 42 PetscFunctionBeginHot; 43 ((DM_Plex *)dm->data)->activePoint = point; 44 PetscFunctionReturn(PETSC_SUCCESS); 45 } 46 47 /* 48 DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point 49 50 Input Parameters: 51 + dm - The output `DM` 52 . ds - The output `DS` 53 . dmIn - The input `DM` 54 . dsIn - The input `DS` 55 . time - The time for this evaluation 56 . fegeom - The FE geometry for this point 57 . fvgeom - The FV geometry for this point 58 . isFE - Flag indicating whether each output field has an FE discretization 59 . sp - The output `PetscDualSpace` for each field 60 . funcs - The evaluation function for each field 61 - ctxs - The user context for each field 62 63 Output Parameter: 64 . values - The value for each dual basis vector in the output dual space 65 66 Level: developer 67 68 .seealso:[](ch_unstructured), `DM`, `DMPLEX`, `PetscDS`, `PetscFEGeom`, `PetscFVCellGeom`, `PetscDualSpace` 69 */ 70 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, PetscScalar values[]) 71 { 72 PetscInt debug = ((DM_Plex *)dm->data)->printProject; 73 PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp; 74 PetscBool isAffine, isCohesive, transform; 75 76 PetscFunctionBeginHot; 77 PetscCall(DMGetCoordinateDim(dmIn, &coordDim)); 78 PetscCall(DMHasBasisTransform(dmIn, &transform)); 79 PetscCall(PetscDSGetNumFields(ds, &Nf)); 80 PetscCall(PetscDSGetComponents(ds, &Nc)); 81 PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 82 /* Get values for closure */ 83 isAffine = fegeom->isAffine; 84 for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 85 void *const ctx = ctxs ? ctxs[f] : NULL; 86 PetscBool cohesive; 87 88 if (!sp[f]) continue; 89 PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 90 PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 91 if (funcs[f]) { 92 if (isFE[f]) { 93 PetscQuadrature allPoints; 94 PetscInt q, dim, numPoints; 95 const PetscReal *points; 96 PetscScalar *pointEval; 97 PetscReal *x; 98 DM rdm; 99 100 PetscCall(PetscDualSpaceGetDM(sp[f], &rdm)); 101 PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 102 PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 103 PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 104 PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x)); 105 PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f])); 106 for (q = 0; q < numPoints; q++, tp++) { 107 const PetscReal *v0; 108 109 if (isAffine) { 110 const PetscReal *refpoint = &points[q * dim]; 111 PetscReal injpoint[3] = {0., 0., 0.}; 112 113 if (dim != fegeom->dim) { 114 if (isCohesive) { 115 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */ 116 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d]; 117 refpoint = injpoint; 118 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim); 119 } 120 CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x); 121 v0 = x; 122 } else { 123 v0 = &fegeom->v[tp * coordDim]; 124 } 125 if (transform) { 126 PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx)); 127 v0 = x; 128 } 129 if (debug > 3) { 130 PetscInt ap; 131 PetscCall(DMPlexGetActivePoint(dm, &ap)); 132 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Project point %" PetscInt_FMT ", analytic: ref (", ap)); 133 for (PetscInt d = 0; d < dim; ++d) { 134 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 135 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)points[q * dim + d])); 136 } 137 PetscCall(PetscPrintf(PETSC_COMM_SELF, ") real (")); 138 for (PetscInt d = 0; d < dim; ++d) { 139 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 140 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)v0[d])); 141 } 142 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 143 } 144 PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx)); 145 } 146 /* Transform point evaluations pointEval[q,c] */ 147 PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval)); 148 PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 149 PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x)); 150 PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 151 v += spDim; 152 if (isCohesive && !cohesive) { 153 for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 154 } 155 } else { 156 for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v])); 157 } 158 } else { 159 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 160 if (isCohesive && !cohesive) { 161 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 162 } 163 } 164 } 165 PetscFunctionReturn(PETSC_SUCCESS); 166 } 167 168 /* 169 DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point 170 171 Input Parameters: 172 + dm - The output DM 173 . ds - The output DS 174 . dmIn - The input DM 175 . dsIn - The input DS 176 . dmAux - The auxiliary DM, which is always for the input space 177 . dsAux - The auxiliary DS, which is always for the input space 178 . time - The time for this evaluation 179 . localU - The local solution 180 . localA - The local auziliary fields 181 . cgeom - The FE geometry for this point 182 . sp - The output PetscDualSpace for each field 183 . p - The point in the output DM 184 . T - Input basis and derivatives for each field tabulated on the quadrature points 185 . TAux - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points 186 . funcs - The evaluation function for each field 187 - ctxs - The user context for each field 188 189 Output Parameter: 190 . values - The value for each dual basis vector in the output dual space 191 192 Level: developer 193 194 Note: 195 Not supported for FV 196 197 .seealso: `DMProjectPoint_Field_Private()` 198 */ 199 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[]) 200 { 201 PetscSection section, sectionAux = NULL; 202 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 203 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 204 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 205 const PetscScalar *constants; 206 PetscReal *x; 207 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2]; 208 PetscFEGeom fegeom, fgeomN[2]; 209 const PetscInt dE = cgeom->dimEmbed, *cone, *ornt; 210 PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 211 PetscBool isAffine, isCohesive, isCohesiveIn, transform; 212 DMPolytopeType qct; 213 214 PetscFunctionBeginHot; 215 PetscCall(PetscDSGetNumFields(ds, &Nf)); 216 PetscCall(PetscDSGetComponents(ds, &Nc)); 217 PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 218 PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 219 PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn)); 220 PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff)); 221 PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x)); 222 PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x)); 223 PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL)); 224 PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants)); 225 PetscCall(DMHasBasisTransform(dmIn, &transform)); 226 PetscCall(DMGetLocalSection(dmIn, §ion)); 227 PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp)); 228 // Get cohesive cell hanging off face 229 if (isCohesiveIn) { 230 PetscCall(DMPlexGetCellType(dmIn, inp, &qct)); 231 if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) { 232 DMPolytopeType ct; 233 const PetscInt *support; 234 PetscInt Ns, s; 235 236 PetscCall(DMPlexGetSupport(dmIn, inp, &support)); 237 PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns)); 238 for (s = 0; s < Ns; ++s) { 239 PetscCall(DMPlexGetCellType(dmIn, support[s], &ct)); 240 if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) { 241 inp = support[s]; 242 break; 243 } 244 } 245 PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp); 246 PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff)); 247 PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt)); 248 face[0] = 0; 249 face[1] = 0; 250 } 251 } 252 if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients)); 253 if (dmAux) { 254 PetscInt subp; 255 256 PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 257 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 258 PetscCall(DMGetLocalSection(dmAux, §ionAux)); 259 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 260 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 261 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 262 PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 263 } 264 /* Get values for closure */ 265 isAffine = cgeom->isAffine; 266 fegeom.dim = cgeom->dim; 267 fegeom.dimEmbed = cgeom->dimEmbed; 268 if (isCohesiveIn) { 269 fgeomN[0].dim = cgeom->dim; 270 fgeomN[0].dimEmbed = cgeom->dimEmbed; 271 fgeomN[1].dim = cgeom->dim; 272 fgeomN[1].dimEmbed = cgeom->dimEmbed; 273 } 274 if (isAffine) { 275 fegeom.v = x; 276 fegeom.xi = cgeom->xi; 277 fegeom.J = cgeom->J; 278 fegeom.invJ = cgeom->invJ; 279 fegeom.detJ = cgeom->detJ; 280 if (isCohesiveIn) { 281 fgeomN[0].J = cgeom->suppJ[0]; 282 fgeomN[0].invJ = cgeom->suppInvJ[0]; 283 fgeomN[0].detJ = cgeom->suppDetJ[0]; 284 fgeomN[1].J = cgeom->suppJ[1]; 285 fgeomN[1].invJ = cgeom->suppInvJ[1]; 286 fgeomN[1].detJ = cgeom->suppDetJ[1]; 287 } 288 } 289 for (f = 0, v = 0; f < Nf; ++f) { 290 PetscQuadrature allPoints; 291 PetscInt q, dim, numPoints; 292 const PetscReal *points; 293 PetscScalar *pointEval; 294 PetscBool cohesive; 295 DM dm; 296 297 if (!sp[f]) continue; 298 PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 299 PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 300 if (!funcs[f]) { 301 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 302 if (isCohesive && !cohesive) { 303 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 304 } 305 continue; 306 } 307 const PetscInt ***perms; 308 PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 309 PetscCall(PetscDualSpaceGetSymmetries(sp[f], &perms, NULL)); 310 PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 311 PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 312 PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 313 for (q = 0; q < numPoints; ++q, ++tp) { 314 PetscInt qpt[2]; 315 316 if (isCohesiveIn) { 317 qpt[0] = perms ? perms[0][ornt[0]][q] : q; 318 qpt[1] = perms ? perms[0][DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0)][q] : q; 319 } 320 if (isAffine) { 321 CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x); 322 } else { 323 fegeom.v = &cgeom->v[tp * dE]; 324 fegeom.J = &cgeom->J[tp * dE * dE]; 325 fegeom.invJ = &cgeom->invJ[tp * dE * dE]; 326 fegeom.detJ = &cgeom->detJ[tp]; 327 if (isCohesiveIn) { 328 fgeomN[0].J = &cgeom->suppJ[0][tp * dE * dE]; 329 fgeomN[0].invJ = &cgeom->suppInvJ[0][tp * dE * dE]; 330 fgeomN[0].detJ = &cgeom->suppDetJ[0][tp]; 331 fgeomN[1].J = &cgeom->suppJ[1][tp * dE * dE]; 332 fgeomN[1].invJ = &cgeom->suppInvJ[1][tp * dE * dE]; 333 fgeomN[1].detJ = &cgeom->suppDetJ[1][tp]; 334 } 335 } 336 if (coefficients) { 337 if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t)); 338 else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t)); 339 } 340 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 341 if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx)); 342 (*funcs[f])(dE, NfIn, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f] * q]); 343 } 344 PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 345 PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 346 v += spDim; 347 /* TODO: For now, set both sides equal, but this should use info from other support cell */ 348 if (isCohesive && !cohesive) { 349 for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 350 } 351 } 352 if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients)); 353 if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 354 if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt)); 355 PetscFunctionReturn(PETSC_SUCCESS); 356 } 357 358 static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[]) 359 { 360 PetscSection section, sectionAux = NULL; 361 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 362 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 363 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 364 const PetscScalar *constants; 365 PetscReal *x; 366 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2]; 367 PetscFEGeom fegeom, cgeom, fgeomN[2]; 368 const PetscInt dE = fgeom->dimEmbed, *cone, *ornt; 369 PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 370 PetscBool isAffine, isCohesive, isCohesiveIn, transform; 371 DMPolytopeType qct; 372 373 PetscFunctionBeginHot; 374 PetscCall(PetscDSGetNumFields(ds, &Nf)); 375 PetscCall(PetscDSGetComponents(ds, &Nc)); 376 PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 377 PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 378 PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn)); 379 PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff)); 380 PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x)); 381 PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x)); 382 PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL)); 383 PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants)); 384 PetscCall(DMHasBasisTransform(dmIn, &transform)); 385 PetscCall(DMGetLocalSection(dmIn, §ion)); 386 PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp)); 387 // Get cohesive cell hanging off face 388 if (isCohesiveIn) { 389 PetscCall(DMPlexGetCellType(dmIn, inp, &qct)); 390 if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) { 391 DMPolytopeType ct; 392 const PetscInt *support; 393 PetscInt Ns, s; 394 395 PetscCall(DMPlexGetSupport(dmIn, inp, &support)); 396 PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns)); 397 for (s = 0; s < Ns; ++s) { 398 PetscCall(DMPlexGetCellType(dmIn, support[s], &ct)); 399 if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) { 400 inp = support[s]; 401 break; 402 } 403 } 404 PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp); 405 PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff)); 406 PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt)); 407 face[0] = 0; 408 face[1] = 0; 409 } 410 } 411 if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients)); 412 if (dmAux) { 413 PetscInt subp; 414 415 PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 416 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 417 PetscCall(DMGetLocalSection(dmAux, §ionAux)); 418 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 419 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 420 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 421 PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 422 } 423 /* Get values for closure */ 424 isAffine = fgeom->isAffine; 425 fegeom.dim = fgeom->dim; 426 fegeom.dimEmbed = fgeom->dimEmbed; 427 fegeom.n = NULL; 428 fegeom.J = NULL; 429 fegeom.v = NULL; 430 fegeom.xi = NULL; 431 cgeom.dim = fgeom->dim; 432 cgeom.dimEmbed = fgeom->dimEmbed; 433 if (isCohesiveIn) { 434 fgeomN[0].dim = fgeom->dim; 435 fgeomN[0].dimEmbed = fgeom->dimEmbed; 436 fgeomN[1].dim = fgeom->dim; 437 fgeomN[1].dimEmbed = fgeom->dimEmbed; 438 } 439 if (isAffine) { 440 fegeom.v = x; 441 fegeom.xi = fgeom->xi; 442 fegeom.J = fgeom->J; 443 fegeom.invJ = fgeom->invJ; 444 fegeom.detJ = fgeom->detJ; 445 fegeom.n = fgeom->n; 446 447 cgeom.J = fgeom->suppJ[0]; 448 cgeom.invJ = fgeom->suppInvJ[0]; 449 cgeom.detJ = fgeom->suppDetJ[0]; 450 451 if (isCohesiveIn) { 452 fgeomN[0].J = fgeom->suppJ[0]; 453 fgeomN[0].invJ = fgeom->suppInvJ[0]; 454 fgeomN[0].detJ = fgeom->suppDetJ[0]; 455 fgeomN[1].J = fgeom->suppJ[1]; 456 fgeomN[1].invJ = fgeom->suppInvJ[1]; 457 fgeomN[1].detJ = fgeom->suppDetJ[1]; 458 } 459 } 460 for (f = 0, v = 0; f < Nf; ++f) { 461 PetscQuadrature allPoints; 462 PetscInt q, dim, numPoints; 463 const PetscReal *points; 464 PetscScalar *pointEval; 465 PetscBool cohesive; 466 DM dm; 467 468 if (!sp[f]) continue; 469 PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 470 PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 471 if (!funcs[f]) { 472 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 473 if (isCohesive && !cohesive) { 474 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 475 } 476 continue; 477 } 478 PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 479 PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 480 PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 481 PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 482 for (q = 0; q < numPoints; ++q, ++tp) { 483 PetscInt qpt[2]; 484 485 if (isCohesiveIn) { 486 // These points are not integration quadratures, but dual space quadratures 487 // If they had multiple points we should match them from both sides, similar to hybrid residual eval 488 qpt[0] = qpt[1] = q; 489 } 490 if (isAffine) { 491 CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x); 492 } else { 493 fegeom.v = &fgeom->v[tp * dE]; 494 fegeom.J = &fgeom->J[tp * dE * dE]; 495 fegeom.invJ = &fgeom->invJ[tp * dE * dE]; 496 fegeom.detJ = &fgeom->detJ[tp]; 497 fegeom.n = &fgeom->n[tp * dE]; 498 499 cgeom.J = &fgeom->suppJ[0][tp * dE * dE]; 500 cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE]; 501 cgeom.detJ = &fgeom->suppDetJ[0][tp]; 502 if (isCohesiveIn) { 503 fgeomN[0].J = &fgeom->suppJ[0][tp * dE * dE]; 504 fgeomN[0].invJ = &fgeom->suppInvJ[0][tp * dE * dE]; 505 fgeomN[0].detJ = &fgeom->suppDetJ[0][tp]; 506 fgeomN[1].J = &fgeom->suppJ[1][tp * dE * dE]; 507 fgeomN[1].invJ = &fgeom->suppInvJ[1][tp * dE * dE]; 508 fgeomN[1].detJ = &fgeom->suppDetJ[1][tp]; 509 } 510 } 511 /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */ 512 if (coefficients) { 513 if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t)); 514 else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t)); 515 } 516 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 517 if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx)); 518 (*funcs[f])(dE, NfIn, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f] * q]); 519 } 520 PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 521 PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 522 v += spDim; 523 /* TODO: For now, set both sides equal, but this should use info from other support cell */ 524 if (isCohesive && !cohesive) { 525 for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 526 } 527 } 528 if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients)); 529 if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 530 if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt)); 531 PetscFunctionReturn(PETSC_SUCCESS); 532 } 533 534 static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscFEGeom *fegeom, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 535 { 536 PetscFVCellGeom fvgeom; 537 PetscInt dim, dimEmbed; 538 539 PetscFunctionBeginHot; 540 PetscCall(DMGetDimension(dm, &dim)); 541 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 542 if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL)); 543 switch (type) { 544 case DM_BC_ESSENTIAL: 545 case DM_BC_NATURAL: 546 PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values)); 547 break; 548 case DM_BC_ESSENTIAL_FIELD: 549 case DM_BC_NATURAL_FIELD: 550 PetscCall(DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values)); 551 break; 552 case DM_BC_ESSENTIAL_BD_FIELD: 553 PetscCall(DMProjectPoint_BdField_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values)); 554 break; 555 default: 556 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type); 557 } 558 PetscFunctionReturn(PETSC_SUCCESS); 559 } 560 561 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 562 { 563 PetscReal *points; 564 PetscInt f, numPoints; 565 566 PetscFunctionBegin; 567 if (!dim) { 568 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 569 PetscFunctionReturn(PETSC_SUCCESS); 570 } 571 numPoints = 0; 572 for (f = 0; f < Nf; ++f) { 573 if (funcs[f]) { 574 PetscQuadrature fAllPoints; 575 PetscInt fNumPoints; 576 577 PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 578 PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL)); 579 numPoints += fNumPoints; 580 } 581 } 582 PetscCall(PetscMalloc1(dim * numPoints, &points)); 583 numPoints = 0; 584 for (f = 0; f < Nf; ++f) { 585 if (funcs[f]) { 586 PetscQuadrature fAllPoints; 587 PetscInt qdim, fNumPoints, q; 588 const PetscReal *fPoints; 589 590 PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 591 PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL)); 592 PetscCheck(qdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %" PetscInt_FMT " for dual basis does not match input dimension %" PetscInt_FMT, qdim, dim); 593 for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q]; 594 numPoints += fNumPoints; 595 } 596 } 597 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 598 PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL)); 599 PetscFunctionReturn(PETSC_SUCCESS); 600 } 601 602 /*@C 603 DMGetFirstLabeledPoint - Find first labeled `point` in `odm` such that the corresponding point in `dm` has the specified `height`. Return `point` and the corresponding `ds`. 604 605 Input Parameters: 606 + dm - the `DM` 607 . odm - the enclosing `DM` 608 . label - label for `DM` domain, or `NULL` for whole domain 609 . numIds - the number of `ids` 610 . ids - An array of the label ids in sequence for the domain 611 - height - Height of target cells in `DMPLEX` topology 612 613 Output Parameters: 614 + point - the first labeled point 615 - ds - the `PetscDS` corresponding to the first labeled point 616 617 Level: developer 618 619 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS` 620 @*/ 621 PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds) 622 { 623 DM plex; 624 DMEnclosureType enc; 625 PetscInt ls = -1; 626 627 PetscFunctionBegin; 628 if (point) *point = -1; 629 if (!label) PetscFunctionReturn(PETSC_SUCCESS); 630 PetscCall(DMGetEnclosureRelation(dm, odm, &enc)); 631 PetscCall(DMConvert(dm, DMPLEX, &plex)); 632 for (PetscInt i = 0; i < numIds; ++i) { 633 IS labelIS; 634 PetscInt num_points, pStart, pEnd; 635 PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS)); 636 if (!labelIS) continue; /* No points with that id on this process */ 637 PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd)); 638 PetscCall(ISGetSize(labelIS, &num_points)); 639 if (num_points) { 640 const PetscInt *points; 641 PetscCall(ISGetIndices(labelIS, &points)); 642 for (PetscInt i = 0; i < num_points; i++) { 643 PetscInt point; 644 PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point)); 645 if (pStart <= point && point < pEnd) { 646 ls = point; 647 if (ds) { 648 // If this is a face of a cohesive cell, then prefer that DS 649 if (height == 1) { 650 const PetscInt *supp; 651 PetscInt suppSize; 652 DMPolytopeType ct; 653 654 PetscCall(DMPlexGetSupport(dm, ls, &supp)); 655 PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize)); 656 for (PetscInt s = 0; s < suppSize; ++s) { 657 PetscCall(DMPlexGetCellType(dm, supp[s], &ct)); 658 if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) { 659 ls = supp[s]; 660 break; 661 } 662 } 663 } 664 PetscCall(DMGetCellDS(dm, ls, ds, NULL)); 665 } 666 if (ls >= 0) break; 667 } 668 } 669 PetscCall(ISRestoreIndices(labelIS, &points)); 670 } 671 PetscCall(ISDestroy(&labelIS)); 672 if (ls >= 0) break; 673 } 674 if (point) *point = ls; 675 PetscCall(DMDestroy(&plex)); 676 PetscFunctionReturn(PETSC_SUCCESS); 677 } 678 679 /* 680 This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 681 682 There are several different scenarios: 683 684 1) Volumetric mesh with volumetric auxiliary data 685 686 Here minHeight=0 since we loop over cells. 687 688 2) Boundary mesh with boundary auxiliary data 689 690 Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation. 691 692 3) Volumetric mesh with boundary auxiliary data 693 694 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. 695 696 4) Volumetric input mesh with boundary output mesh 697 698 Here we must get a subspace for the input DS 699 700 The maxHeight is used to support enforcement of constraints in DMForest. 701 702 If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 703 704 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. 705 - 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 706 is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When the effective height is nonzero, we need to extract 707 dual spaces for the boundary from our input spaces. 708 - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 709 710 We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 711 712 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. 713 */ 714 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, InsertMode mode, Vec localX) 715 { 716 DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 717 DMEnclosureType encIn, encAux; 718 PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 719 Vec localA = NULL, tv; 720 IS fieldIS; 721 PetscSection section; 722 PetscDualSpace *sp, *cellsp, *spIn, *cellspIn; 723 PetscTabulation *T = NULL, *TAux = NULL; 724 PetscInt *Nc; 725 PetscInt dim, dimEmbed, depth, pStart, pEnd, lStart = PETSC_DETERMINE, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f; 726 PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform; 727 DMField coordField; 728 DMLabel depthLabel; 729 PetscQuadrature allPoints = NULL; 730 731 PetscFunctionBegin; 732 if (localU) PetscCall(VecGetDM(localU, &dmIn)); 733 else dmIn = dm; 734 PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA)); 735 if (localA) PetscCall(VecGetDM(localA, &dmAux)); 736 else dmAux = NULL; 737 PetscCall(DMConvert(dm, DMPLEX, &plex)); 738 PetscCall(DMConvert(dmIn, DMPLEX, &plexIn)); 739 PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn)); 740 PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux)); 741 PetscCall(DMGetDimension(dm, &dim)); 742 PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight)); 743 PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm)); 744 PetscCall(DMGetBasisTransformVec_Internal(dm, &tv)); 745 PetscCall(DMHasBasisTransform(dm, &transform)); 746 /* Auxiliary information can only be used with interpolation of field functions */ 747 if (dmAux) { 748 PetscCall(DMConvert(dmAux, DMPLEX, &plexAux)); 749 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) PetscCheck(localA, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector"); 750 } 751 if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL)); 752 PetscCall(DMGetCoordinateField(dm, &coordField)); 753 PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field"); 754 /**** No collective calls below this point ****/ 755 /* Determine height for iteration of all meshes */ 756 { 757 DMPolytopeType ct, ctIn, ctAux; 758 PetscInt p, pStartIn, pStartAux, pEndAux; 759 PetscInt dim = -1, dimIn = -1, dimAux = -1; 760 761 PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd)); 762 if (pEnd > pStart) { 763 PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL)); 764 p = lStart < 0 ? pStart : lStart; 765 PetscCall(DMPlexGetCellType(plex, p, &ct)); 766 dim = DMPolytopeTypeGetDim(ct); 767 PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn)); 768 PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL)); 769 PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn)); 770 dimIn = DMPolytopeTypeGetDim(ctIn); 771 if (dmAux) { 772 PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux)); 773 PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux)); 774 if (pStartAux < pEndAux) { 775 PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux)); 776 dimAux = DMPolytopeTypeGetDim(ctAux); 777 } 778 } else dimAux = dim; 779 } else { 780 PetscCall(DMDestroy(&plex)); 781 PetscCall(DMDestroy(&plexIn)); 782 if (dmAux) PetscCall(DMDestroy(&plexAux)); 783 PetscFunctionReturn(PETSC_SUCCESS); 784 } 785 if (dim < 0) { 786 DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL; 787 788 /* Fall back to determination based on being a submesh */ 789 PetscCall(DMPlexGetSubpointMap(plex, &spmap)); 790 PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn)); 791 if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux)); 792 dim = spmap ? 1 : 0; 793 dimIn = spmapIn ? 1 : 0; 794 dimAux = spmapAux ? 1 : 0; 795 } 796 { 797 PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), dimAux < 0 ? PETSC_INT_MAX : dimAux); 798 PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux; 799 800 PetscCheck(PetscAbsInt(dimProj - dim) <= 1 && PetscAbsInt(dimProj - dimIn) <= 1 && PetscAbsInt(dimProj - dimAuxEff) <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension"); 801 if (dimProj < dim) minHeight = 1; 802 htInc = dim - dimProj; 803 htIncIn = dimIn - dimProj; 804 htIncAux = dimAuxEff - dimProj; 805 } 806 } 807 PetscCall(DMPlexGetDepth(plex, &depth)); 808 PetscCall(DMPlexGetDepthLabel(plex, &depthLabel)); 809 PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight)); 810 maxHeight = PetscMax(maxHeight, minHeight); 811 PetscCheck(maxHeight >= 0 && maxHeight <= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", maxHeight, dim); 812 PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds)); 813 if (!ds) PetscCall(DMGetDS(dm, &ds)); 814 PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn)); 815 if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn)); 816 if (!dsIn) { 817 if (encIn == DM_ENC_SUPERMESH) { 818 PetscInt p = pStart, pIn; 819 820 PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn)); 821 // If the input mesh is higher dimensional than the output mesh, get a cell from the output mesh 822 if (htIncIn) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL)); 823 PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn)); 824 PetscCall(DMGetCellDS(dmIn, pIn, &dsIn, NULL)); 825 } else PetscCall(DMGetDS(dmIn, &dsIn)); 826 } 827 PetscCall(PetscDSGetNumFields(ds, &Nf)); 828 PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 829 PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn)); 830 if (isCohesiveIn) --htIncIn; // Should be rearranged 831 PetscCall(DMGetNumFields(dm, &NfTot)); 832 PetscCall(DMFindRegionNum(dm, ds, ®ionNum)); 833 PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL)); 834 PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 835 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 836 PetscCall(DMGetLocalSection(dm, §ion)); 837 if (dmAux) { 838 if (encAux == DM_ENC_SUPERMESH) { 839 PetscInt p = pStart, pAux; 840 841 // If the auxiliary mesh is higher dimensional than the output mesh, get a cell from the output mesh 842 if (htIncAux) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL)); 843 PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, lStart < 0 ? p : lStart, &pAux)); 844 PetscCall(DMGetCellDS(dmAux, pAux, &dsAux, NULL)); 845 } else PetscCall(DMGetDS(dmAux, &dsAux)); 846 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 847 } 848 PetscCall(PetscDSGetComponents(ds, &Nc)); 849 PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn)); 850 if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn)); 851 else { 852 cellsp = sp; 853 cellspIn = spIn; 854 } 855 /* Get cell dual spaces */ 856 for (f = 0; f < Nf; ++f) { 857 PetscDiscType disctype; 858 859 PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype)); 860 if (disctype == PETSC_DISC_FE) { 861 PetscFE fe; 862 863 isFE[f] = PETSC_TRUE; 864 hasFE = PETSC_TRUE; 865 PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 866 PetscCall(PetscFEGetDualSpace(fe, &cellsp[f])); 867 } else if (disctype == PETSC_DISC_FV) { 868 PetscFV fv; 869 870 isFE[f] = PETSC_FALSE; 871 hasFV = PETSC_TRUE; 872 PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv)); 873 PetscCall(PetscFVGetDualSpace(fv, &cellsp[f])); 874 } else { 875 isFE[f] = PETSC_FALSE; 876 cellsp[f] = NULL; 877 } 878 } 879 for (f = 0; f < NfIn; ++f) { 880 PetscDiscType disctype; 881 882 PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 883 if (disctype == PETSC_DISC_FE) { 884 PetscFE fe; 885 886 PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe)); 887 PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f])); 888 } else if (disctype == PETSC_DISC_FV) { 889 PetscFV fv; 890 891 PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv)); 892 PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f])); 893 } else { 894 cellspIn[f] = NULL; 895 } 896 } 897 for (f = 0; f < Nf; ++f) { 898 if (!htInc) { 899 sp[f] = cellsp[f]; 900 } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f])); 901 } 902 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 903 PetscFE fem, subfem; 904 PetscDiscType disctype; 905 const PetscReal *points; 906 PetscInt numPoints, k; 907 908 PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 909 PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints)); 910 PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL)); 911 PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux)); 912 for (f = 0; f < NfIn; ++f) { 913 if (!htIncIn) { 914 spIn[f] = cellspIn[f]; 915 } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f])); 916 917 PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 918 if (disctype != PETSC_DISC_FE) continue; 919 PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem)); 920 if (!htIncIn) { 921 subfem = fem; 922 } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem)); 923 PetscCall(PetscDSGetJetDegree(dsIn, f, &k)); 924 PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &T[f])); 925 } 926 for (f = 0; f < NfAux; ++f) { 927 PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype)); 928 if (disctype != PETSC_DISC_FE) continue; 929 PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem)); 930 if (!htIncAux) { 931 subfem = fem; 932 } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem)); 933 PetscCall(PetscDSGetJetDegree(dsAux, f, &k)); 934 PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &TAux[f])); 935 } 936 } 937 /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 938 for (h = minHeight; h <= maxHeight; h++) { 939 PetscInt hEff = h - minHeight + htInc; 940 PetscInt hEffIn = h - minHeight + htIncIn; 941 PetscInt hEffAux = h - minHeight + htIncAux; 942 PetscDS dsEff = ds; 943 PetscDS dsEffIn = dsIn; 944 PetscDS dsEffAux = dsAux; 945 PetscScalar *values; 946 PetscBool *fieldActive; 947 PetscInt maxDegree; 948 PetscInt p, spDim, totDim, numValues; 949 IS heightIS; 950 951 if (h > minHeight) { 952 for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f])); 953 } 954 PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd)); 955 PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL)); 956 PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS)); 957 if (pEnd <= pStart) { 958 PetscCall(ISDestroy(&heightIS)); 959 continue; 960 } 961 /* Compute totDim, the number of dofs in the closure of a point at this height */ 962 totDim = 0; 963 for (f = 0; f < Nf; ++f) { 964 PetscBool cohesive; 965 966 if (!sp[f]) continue; 967 PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 968 PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 969 totDim += spDim; 970 if (isCohesive && !cohesive) totDim += spDim; 971 } 972 p = lStart < 0 ? pStart : lStart; 973 PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL)); 974 PetscCheck(numValues == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT " in [%" PetscInt_FMT ", %" PetscInt_FMT "]", p, numValues, totDim, h, minHeight, maxHeight); 975 if (!totDim) { 976 PetscCall(ISDestroy(&heightIS)); 977 continue; 978 } 979 if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff)); 980 /* Compute totDimIn, the number of dofs in the closure of a point at this height */ 981 if (localU) { 982 PetscInt totDimIn, pIn, numValuesIn; 983 984 totDimIn = 0; 985 for (f = 0; f < NfIn; ++f) { 986 PetscBool cohesive; 987 988 if (!spIn[f]) continue; 989 PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive)); 990 PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim)); 991 totDimIn += spDim; 992 if (isCohesiveIn && !cohesive) totDimIn += spDim; 993 } 994 PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn)); 995 PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL)); 996 // TODO We could check that pIn is a cohesive cell for this check 997 PetscCheck(isCohesiveIn || (numValuesIn == totDimIn), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT, pIn, numValuesIn, totDimIn, htIncIn); 998 if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn)); 999 } 1000 if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux)); 1001 /* Loop over points at this height */ 1002 PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values)); 1003 PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive)); 1004 { 1005 const PetscInt *fields; 1006 1007 PetscCall(ISGetIndices(fieldIS, &fields)); 1008 for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE; 1009 for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 1010 PetscCall(ISRestoreIndices(fieldIS, &fields)); 1011 } 1012 if (label) { 1013 PetscInt i; 1014 1015 for (i = 0; i < numIds; ++i) { 1016 IS pointIS, isectIS; 1017 const PetscInt *points; 1018 PetscInt n; 1019 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 1020 PetscQuadrature quad = NULL; 1021 1022 PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS)); 1023 if (!pointIS) continue; /* No points with that id on this process */ 1024 PetscCall(ISIntersect(pointIS, heightIS, &isectIS)); 1025 PetscCall(ISDestroy(&pointIS)); 1026 if (!isectIS) continue; 1027 PetscCall(ISGetLocalSize(isectIS, &n)); 1028 PetscCall(ISGetIndices(isectIS, &points)); 1029 PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree)); 1030 if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad)); 1031 if (!quad) { 1032 if (!h && allPoints) { 1033 quad = allPoints; 1034 allPoints = NULL; 1035 } else { 1036 PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad)); 1037 } 1038 } 1039 PetscFEGeomMode geommode = htInc && h == minHeight ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC; 1040 1041 if (n) { 1042 PetscInt depth, dep; 1043 1044 PetscCall(DMPlexGetDepth(dm, &depth)); 1045 PetscCall(DMPlexGetPointDepth(dm, points[0], &dep)); 1046 if (dep < depth && h == minHeight) geommode = PETSC_FEGEOM_BOUNDARY; 1047 } 1048 PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, geommode, &fegeom)); 1049 for (p = 0; p < n; ++p) { 1050 const PetscInt point = points[p]; 1051 1052 PetscCall(PetscArrayzero(values, numValues)); 1053 PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom)); 1054 PetscCall(DMPlexSetActivePoint(dm, point)); 1055 PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values)); 1056 if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values)); 1057 PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode)); 1058 } 1059 PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom)); 1060 PetscCall(PetscFEGeomDestroy(&fegeom)); 1061 PetscCall(PetscQuadratureDestroy(&quad)); 1062 PetscCall(ISRestoreIndices(isectIS, &points)); 1063 PetscCall(ISDestroy(&isectIS)); 1064 } 1065 } else { 1066 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 1067 PetscQuadrature quad = NULL; 1068 IS pointIS; 1069 1070 PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS)); 1071 PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree)); 1072 if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad)); 1073 if (!quad) { 1074 if (!h && allPoints) { 1075 quad = allPoints; 1076 allPoints = NULL; 1077 } else { 1078 PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad)); 1079 } 1080 } 1081 PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC, &fegeom)); 1082 for (p = pStart; p < pEnd; ++p) { 1083 PetscCall(PetscArrayzero(values, numValues)); 1084 PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom)); 1085 PetscCall(DMPlexSetActivePoint(dm, p)); 1086 PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values)); 1087 if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values)); 1088 PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode)); 1089 } 1090 PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom)); 1091 PetscCall(PetscFEGeomDestroy(&fegeom)); 1092 PetscCall(PetscQuadratureDestroy(&quad)); 1093 PetscCall(ISDestroy(&pointIS)); 1094 } 1095 PetscCall(ISDestroy(&heightIS)); 1096 PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values)); 1097 PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive)); 1098 } 1099 /* Cleanup */ 1100 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 1101 for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f])); 1102 for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f])); 1103 PetscCall(PetscFree2(T, TAux)); 1104 } 1105 PetscCall(PetscQuadratureDestroy(&allPoints)); 1106 PetscCall(PetscFree3(isFE, sp, spIn)); 1107 if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn)); 1108 PetscCall(DMDestroy(&plex)); 1109 PetscCall(DMDestroy(&plexIn)); 1110 if (dmAux) PetscCall(DMDestroy(&plexAux)); 1111 PetscFunctionReturn(PETSC_SUCCESS); 1112 } 1113 1114 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 1115 { 1116 PetscFunctionBegin; 1117 PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 1118 PetscFunctionReturn(PETSC_SUCCESS); 1119 } 1120 1121 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) 1122 { 1123 PetscFunctionBegin; 1124 PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 1125 PetscFunctionReturn(PETSC_SUCCESS); 1126 } 1127 1128 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX) 1129 { 1130 PetscFunctionBegin; 1131 PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 1132 PetscFunctionReturn(PETSC_SUCCESS); 1133 } 1134 1135 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX) 1136 { 1137 PetscFunctionBegin; 1138 PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 1139 PetscFunctionReturn(PETSC_SUCCESS); 1140 } 1141 1142 PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX) 1143 { 1144 PetscFunctionBegin; 1145 PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 1146 PetscFunctionReturn(PETSC_SUCCESS); 1147 } 1148