1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petsc/private/petscfeimpl.h> 4 5 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS prob, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], 6 PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 7 PetscScalar values[]) 8 { 9 PetscInt coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp; 10 PetscBool isAffine; 11 PetscErrorCode ierr; 12 13 PetscFunctionBeginHot; 14 ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr); 15 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 16 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 17 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 18 /* Get values for closure */ 19 isAffine = fegeom->isAffine; 20 for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 21 void * const ctx = ctxs ? ctxs[f] : NULL; 22 23 if (!sp[f]) continue; 24 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 25 if (funcs[f]) { 26 if (isFE[f]) { 27 PetscQuadrature allPoints; 28 PetscInt q, dim, numPoints; 29 const PetscReal *points; 30 PetscScalar *pointEval; 31 PetscReal *x; 32 DM dm; 33 34 ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 35 ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 36 ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 37 ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 38 ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 39 ierr = DMGetWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 40 for (q = 0; q < numPoints; q++, tp++) { 41 const PetscReal *v0; 42 43 if (isAffine) { 44 CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x); 45 v0 = x; 46 } else { 47 v0 = &fegeom->v[tp*coordDim]; 48 } 49 ierr = (*funcs[f])(coordDim,time,v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 50 } 51 ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 52 ierr = DMRestoreWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 53 ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 54 v += spDim; 55 } else { 56 for (d = 0; d < spDim; ++d, ++v) { 57 ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 58 } 59 } 60 } else { 61 for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;} 62 } 63 } 64 PetscFunctionReturn(0); 65 } 66 67 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 68 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 69 void (**funcs)(PetscInt, PetscInt, PetscInt, 70 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 71 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 72 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 73 PetscScalar values[]) 74 { 75 PetscSection section, sectionAux = NULL; 76 PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc; 77 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 78 PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 79 const PetscScalar *constants; 80 PetscReal *x; 81 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 82 const PetscInt dE = fegeom->dimEmbed; 83 PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 84 PetscBool isAffine; 85 PetscErrorCode ierr; 86 87 PetscFunctionBeginHot; 88 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 89 ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 90 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 91 ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 92 ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 93 ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 94 ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 95 ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 96 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 97 ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 98 if (dmAux) { 99 PetscInt subp; 100 101 ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr); 102 ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 103 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 104 ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 105 ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 106 ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 107 ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 108 ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 109 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 110 ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 111 ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 112 } 113 /* Get values for closure */ 114 isAffine = fegeom->isAffine; 115 for (f = 0, v = 0; f < Nf; ++f) { 116 PetscQuadrature allPoints; 117 PetscInt q, dim, numPoints; 118 const PetscReal *points; 119 PetscScalar *pointEval; 120 DM dm; 121 122 if (!sp[f]) continue; 123 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 124 if (!funcs[f]) { 125 for (d = 0; d < spDim; d++, v++) values[v] = 0.; 126 continue; 127 } 128 ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 129 ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 130 ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 131 ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 132 for (q = 0; q < numPoints; ++q, ++tp) { 133 const PetscReal *v0; 134 const PetscReal *invJ; 135 136 if (isAffine) { 137 CoordinatesRefToReal(dE, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x); 138 v0 = x; 139 invJ = fegeom->invJ; 140 } else { 141 v0 = &fegeom->v[tp*dE]; 142 invJ = &fegeom->invJ[tp*dE*dE]; 143 } 144 EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t); 145 if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t); 146 (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, v0, numConstants, constants, &pointEval[Nc[f]*q]); 147 } 148 ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 149 ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 150 v += spDim; 151 } 152 ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 153 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 154 PetscFunctionReturn(0); 155 } 156 157 static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, PetscFEGeom *fegeom, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], 158 PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 159 PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 160 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 161 { 162 PetscFVCellGeom fvgeom; 163 PetscInt dim, dimEmbed; 164 PetscErrorCode ierr; 165 166 PetscFunctionBeginHot; 167 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 168 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 169 if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 170 switch (type) { 171 case DM_BC_ESSENTIAL: 172 case DM_BC_NATURAL: 173 ierr = DMProjectPoint_Func_Private(dm, prob, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break; 174 case DM_BC_ESSENTIAL_FIELD: 175 case DM_BC_NATURAL_FIELD: 176 ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps, 177 basisTab, basisDerTab, basisTabAux, basisDerTabAux, 178 (void (**)(PetscInt, PetscInt, PetscInt, 179 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 180 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 181 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 182 default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 183 } 184 PetscFunctionReturn(0); 185 } 186 187 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 188 { 189 PetscReal *points; 190 PetscInt f, numPoints; 191 PetscErrorCode ierr; 192 193 PetscFunctionBegin; 194 numPoints = 0; 195 for (f = 0; f < Nf; ++f) { 196 if (funcs[f]) { 197 PetscQuadrature fAllPoints; 198 PetscInt fNumPoints; 199 200 ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 201 ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 202 numPoints += fNumPoints; 203 } 204 } 205 ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 206 numPoints = 0; 207 for (f = 0; f < Nf; ++f) { 208 PetscInt spDim; 209 210 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 211 if (funcs[f]) { 212 PetscQuadrature fAllPoints; 213 PetscInt qdim, fNumPoints, q; 214 const PetscReal *fPoints; 215 216 ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 217 ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 218 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); 219 for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 220 numPoints += fNumPoints; 221 } 222 } 223 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 224 ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 225 PetscFunctionReturn(0); 226 } 227 228 static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob) 229 { 230 DMLabel depthLabel; 231 PetscInt dim, cdepth, ls = -1, i; 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 if (lStart) *lStart = -1; 236 if (!label) PetscFunctionReturn(0); 237 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 238 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 239 cdepth = dim - height; 240 for (i = 0; i < numIds; ++i) { 241 IS pointIS; 242 const PetscInt *points; 243 PetscInt pdepth; 244 245 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 246 if (!pointIS) continue; /* No points with that id on this process */ 247 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 248 ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr); 249 if (pdepth == cdepth) { 250 ls = points[0]; 251 if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);} 252 } 253 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 254 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 255 if (ls >= 0) break; 256 } 257 if (lStart) *lStart = ls; 258 PetscFunctionReturn(0); 259 } 260 261 /* 262 This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 263 264 There are several different scenarios: 265 266 1) Volumetric mesh with volumetric auxiliary data 267 268 Here minHeight=0 since we loop over cells. 269 270 2) Boundary mesh with boundary auxiliary data 271 272 Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation. 273 274 3) Volumetric mesh with boundary auxiliary data 275 276 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. 277 278 The maxHeight is used to support enforcement of constraints in DMForest. 279 280 If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 281 282 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. 283 - 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 284 is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract 285 dual spaces for the boundary from our input spaces. 286 - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 287 288 We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 289 290 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. 291 */ 292 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 293 PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 294 DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 295 InsertMode mode, Vec localX) 296 { 297 DM dmAux = NULL; 298 PetscDS prob = NULL, probAux = NULL; 299 Vec localA = NULL; 300 PetscSection section; 301 PetscDualSpace *sp, *cellsp; 302 PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 303 PetscInt *Nc; 304 PetscInt dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f; 305 PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE; 306 DMField coordField; 307 DMLabel depthLabel; 308 PetscQuadrature allPoints = NULL; 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 313 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 314 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 315 ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr); 316 /* Auxiliary information can only be used with interpolation of field functions */ 317 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 318 if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector"); 319 if (!minHeight && dmAux) { 320 DMLabel spmap; 321 322 /* If dmAux is a surface, then force the projection to take place over a surface */ 323 ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr); 324 if (spmap) { 325 ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr); 326 auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE; 327 } 328 } 329 } 330 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 331 ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr); 332 ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr); 333 maxHeight = PetscMax(maxHeight, minHeight); 334 if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 335 ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr); 336 if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);} 337 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 338 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 339 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 340 if (dmAux) { 341 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 342 ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 343 } 344 ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 345 ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 346 if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 347 else {cellsp = sp;} 348 if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 349 /* Get cell dual spaces */ 350 for (f = 0; f < Nf; ++f) { 351 PetscObject obj; 352 PetscClassId id; 353 354 ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 355 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 356 if (id == PETSCFE_CLASSID) { 357 PetscFE fe = (PetscFE) obj; 358 359 hasFE = PETSC_TRUE; 360 isFE[f] = PETSC_TRUE; 361 ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 362 } else if (id == PETSCFV_CLASSID) { 363 PetscFV fv = (PetscFV) obj; 364 365 hasFV = PETSC_TRUE; 366 isFE[f] = PETSC_FALSE; 367 ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 368 } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 369 } 370 ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 371 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 372 PetscInt effectiveHeight = auxBd ? minHeight : 0; 373 PetscFE fem, subfem; 374 const PetscReal *points; 375 PetscInt numPoints; 376 377 if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 378 for (f = 0; f < Nf; ++f) { 379 if (!effectiveHeight) {sp[f] = cellsp[f];} 380 else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 381 } 382 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr); 383 ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 384 ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 385 for (f = 0; f < Nf; ++f) { 386 if (!isFE[f]) continue; 387 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 388 if (!effectiveHeight) {subfem = fem;} 389 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 390 ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 391 } 392 for (f = 0; f < NfAux; ++f) { 393 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 394 if (!effectiveHeight || auxBd) {subfem = fem;} 395 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 396 ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 397 } 398 } 399 /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 400 for (h = minHeight; h <= maxHeight; h++) { 401 PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 402 PetscDS probEff = prob; 403 PetscScalar *values; 404 PetscBool *fieldActive; 405 PetscInt maxDegree; 406 PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 407 IS heightIS; 408 409 ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 410 ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 411 ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 412 if (!h) { 413 PetscInt cEndInterior; 414 415 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 416 pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 417 } 418 if (pEnd <= pStart) { 419 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 420 continue; 421 } 422 /* Compute totDim, the number of dofs in the closure of a point at this height */ 423 totDim = 0; 424 for (f = 0; f < Nf; ++f) { 425 if (!effectiveHeight) { 426 sp[f] = cellsp[f]; 427 } else { 428 ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 429 if (!sp[f]) continue; 430 } 431 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 432 totDim += spDim; 433 } 434 ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr); 435 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim); 436 if (!totDim) { 437 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 438 continue; 439 } 440 if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);} 441 /* Loop over points at this height */ 442 ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 443 ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 444 for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 445 if (label) { 446 PetscInt i; 447 448 for (i = 0; i < numIds; ++i) { 449 IS pointIS, isectIS; 450 const PetscInt *points; 451 PetscInt n; 452 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 453 PetscQuadrature quad = NULL; 454 455 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 456 if (!pointIS) continue; /* No points with that id on this process */ 457 ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 458 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 459 if (!isectIS) continue; 460 ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 461 ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 462 ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 463 if (maxDegree <= 1) { 464 ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 465 } 466 if (!quad) { 467 if (!h && allPoints) { 468 quad = allPoints; 469 allPoints = NULL; 470 } else { 471 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 472 } 473 } 474 ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 475 for (p = 0; p < n; ++p) { 476 const PetscInt point = points[p]; 477 478 ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 479 ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 480 ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 481 if (ierr) { 482 PetscErrorCode ierr2; 483 ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 484 ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 485 CHKERRQ(ierr); 486 } 487 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr); 488 } 489 ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 490 ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 491 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 492 ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 493 ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 494 } 495 } else { 496 PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 497 PetscQuadrature quad = NULL; 498 IS pointIS; 499 500 ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 501 ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 502 if (maxDegree <= 1) { 503 ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 504 } 505 if (!quad) { 506 if (!h && allPoints) { 507 quad = allPoints; 508 allPoints = NULL; 509 } else { 510 ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 511 } 512 } 513 ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 514 for (p = pStart; p < pEnd; ++p) { 515 ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 516 ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 517 ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 518 if (ierr) { 519 PetscErrorCode ierr2; 520 ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 521 ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 522 CHKERRQ(ierr); 523 } 524 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr); 525 } 526 ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 527 ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 528 ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 529 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 530 } 531 ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 532 ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 533 ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 534 } 535 /* Cleanup */ 536 if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 537 PetscInt effectiveHeight = auxBd ? minHeight : 0; 538 PetscFE fem, subfem; 539 540 for (f = 0; f < Nf; ++f) { 541 if (!isFE[f]) continue; 542 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 543 if (!effectiveHeight) {subfem = fem;} 544 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 545 ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 546 } 547 for (f = 0; f < NfAux; ++f) { 548 ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 549 if (!effectiveHeight || auxBd) {subfem = fem;} 550 else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 551 ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 552 } 553 ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 554 } 555 ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 556 ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 557 if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 558 PetscFunctionReturn(0); 559 } 560 561 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 562 { 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 567 PetscFunctionReturn(0); 568 } 569 570 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) 571 { 572 PetscErrorCode ierr; 573 574 PetscFunctionBegin; 575 ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 580 void (**funcs)(PetscInt, PetscInt, PetscInt, 581 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 582 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 583 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 584 InsertMode mode, Vec localX) 585 { 586 PetscErrorCode ierr; 587 588 PetscFunctionBegin; 589 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 590 PetscFunctionReturn(0); 591 } 592 593 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 594 void (**funcs)(PetscInt, PetscInt, PetscInt, 595 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 596 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 597 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 598 InsertMode mode, Vec localX) 599 { 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606