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