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