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