/*@C DMPlexTSGetGradientDM - Return gradient data layout Input Parameters: + dm - The DM - fv - The PetscFV Output Parameter: . dmGrad - The layout for gradient values Level: developer .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal() @*/ PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad) { DMTS dmts; PetscObject obj; PetscBool computeGradients; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(dm,DM_CLASSID,1); PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2); PetscValidPointer(dmGrad,3); ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr); if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);} ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr); if (!obj) { DM dmGrad; Vec faceGeometry, cellGeometry; ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr); ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr); ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); } ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@C DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function Logically Collective Input Arguments: + dm - DM to associate callback with . riemann - Riemann solver - ctx - optional context for Riemann solve Calling sequence for riemann: $ riemann(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) + x - The coordinates at a point on the interface . n - The normal vector to the interface . uL - The state vector to the left of the interface . uR - The state vector to the right of the interface . flux - output array of flux through the interface - ctx - optional user context Level: beginner .seealso: DMTSSetRHSFunctionLocal() @*/ PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, void (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx), void *ctx) { DMTS dmts; DMTS_Plex *dmplexts; PetscFV fvm; PetscInt Nf; PetscBool computeGradients; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(dm, DM_CLASSID, 1); ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr); ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr); dmplexts->riemann = riemann; dmplexts->rhsfunctionlocalctx = ctx; ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); ierr = DMPlexTSSetupGeometry(dm, fvm, dmplexts);CHKERRQ(ierr); ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmplexts);CHKERRQ(ierr);} ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMPlex, dmplexts);CHKERRQ(ierr); PetscFunctionReturn(0); }
PETSC_EXTERN void PETSC_STDCALL petscfvgetcomputegradients_(PetscFV fvm,PetscBool *computeGradients, int *__ierr ){ *__ierr = PetscFVGetComputeGradients( (PetscFV)PetscToPointer((fvm) ),computeGradients); }
static PetscErrorCode TSComputeRHSFunction_DMPlex(TS ts, PetscReal time, Vec X, Vec F, void *ctx) { DM dm; DMTS_Plex *dmplexts = (DMTS_Plex *) ctx; void (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *) = dmplexts->riemann; PetscFV fvm; PetscLimiter lim; Vec faceGeometry = dmplexts->facegeom; Vec cellGeometry = dmplexts->cellgeom; Vec Grad = NULL, locGrad, locX; DM dmFace, dmCell; DMLabel ghostLabel; PetscCellGeometry fgeom, cgeom; const PetscScalar *facegeom, *cellgeom, *x, *lgrad; PetscScalar *grad, *f, *uL, *uR, *fluxL, *fluxR; PetscReal *centroid, *normal, *vol, *cellPhi; PetscBool computeGradients; PetscInt Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(ts,TS_CLASSID,1); PetscValidHeaderSpecific(X,VEC_CLASSID,3); PetscValidHeaderSpecific(F,VEC_CLASSID,5); ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); ierr = VecZeroEntries(locX);CHKERRQ(ierr); ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); ierr = VecZeroEntries(F);CHKERRQ(ierr); ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); if (computeGradients) { ierr = DMGetGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr); ierr = VecZeroEntries(Grad);CHKERRQ(ierr); ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr); } ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); /* Count faces and reconstruct gradients */ for (face = fStart; face < fEnd; ++face) { const PetscInt *cells; const FaceGeom *fg; const PetscScalar *cx[2]; PetscScalar *cgrad[2]; PetscBool boundary; PetscInt ghost, c, pd, d; ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); if (ghost >= 0) continue; ++numFaces; if (!computeGradients) continue; ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); if (boundary) continue; ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); for (c = 0; c < 2; ++c) { ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr); } for (pd = 0; pd < pdim; ++pd) { PetscScalar delta = cx[1][pd] - cx[0][pd]; for (d = 0; d < dim; ++d) { if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; } } } /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { const PetscInt *faces; const PetscScalar *cx; const CellGeom *cg; PetscScalar *cgrad; PetscInt coneSize, f, pd, d; ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cell, grad, &cgrad);CHKERRQ(ierr); if (!cgrad) continue; /* Unowned overlap cell, we do not compute */ /* Limiter will be minimum value over all neighbors */ for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL; for (f = 0; f < coneSize; ++f) { const PetscScalar *ncx; const CellGeom *ncg; const PetscInt *fcells; PetscInt face = faces[f], ncell, ghost; PetscReal v[3]; PetscBool boundary; ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); if ((ghost >= 0) || boundary) continue; ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); ncell = cell == fcells[0] ? fcells[1] : fcells[0]; ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); WaxpyD(dim, -1, cg->centroid, ncg->centroid, v); for (d = 0; d < pdim; ++d) { /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DotD(dim, &cgrad[d*dim], v); ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); cellPhi[d] = PetscMin(cellPhi[d], phi); } } /* Apply limiter to gradient */ for (pd = 0; pd < pdim; ++pd) /* Scalar limiter applied to each component separately */ for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; } ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad, dmplexts);CHKERRQ(ierr); if (computeGradients) { ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr); ierr = DMGetLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr); ierr = DMGlobalToLocalBegin(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr); ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); } ierr = PetscMalloc7(numFaces*dim,¢roid,numFaces*dim,&normal,numFaces*2,&vol,numFaces*pdim,&uL,numFaces*pdim,&uR,numFaces*pdim,&fluxL,numFaces*pdim,&fluxR);CHKERRQ(ierr); /* Read out values */ for (face = fStart, iface = 0; face < fEnd; ++face) { const PetscInt *cells; const FaceGeom *fg; const CellGeom *cgL, *cgR; const PetscScalar *xL, *xR, *gL, *gR; PetscInt ghost, d; ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); if (ghost >= 0) continue; ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr); if (computeGradients) { PetscReal dxL[3], dxR[3]; ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); WaxpyD(dim, -1, cgL->centroid, fg->centroid, dxL); WaxpyD(dim, -1, cgR->centroid, fg->centroid, dxR); for (d = 0; d < pdim; ++d) { uL[iface*pdim+d] = xL[d] + DotD(dim, &gL[d*dim], dxL); uR[iface*pdim+d] = xR[d] + DotD(dim, &gR[d*dim], dxR); } } else { for (d = 0; d < pdim; ++d) { uL[iface*pdim+d] = xL[d]; uR[iface*pdim+d] = xR[d]; } } for (d = 0; d < dim; ++d) { centroid[iface*dim+d] = fg->centroid[d]; normal[iface*dim+d] = fg->normal[d]; } vol[iface*2+0] = cgL->volume; vol[iface*2+1] = cgR->volume; ++iface; } if (computeGradients) { ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr); ierr = DMRestoreLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr); } ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); fgeom.v0 = centroid; fgeom.n = normal; cgeom.vol = vol; /* Riemann solve */ ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, dmplexts->rhsfunctionlocalctx);CHKERRQ(ierr); /* Insert fluxes */ ierr = VecGetArray(F, &f);CHKERRQ(ierr); for (face = fStart, iface = 0; face < fEnd; ++face) { const PetscInt *cells; PetscScalar *fL, *fR; PetscInt ghost, d; ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); if (ghost >= 0) continue; ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr); ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr); for (d = 0; d < pdim; ++d) { if (fL) fL[d] -= fluxL[iface*pdim+d]; if (fR) fR[d] += fluxR[iface*pdim+d]; } ++iface; } ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr); ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); PetscFunctionReturn(0); }