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