chombo-discharge
|
Godunov implementation for advection. More...
#include <CD_CdrGodunov.H>
Public Member Functions | |
CdrGodunov () | |
Constructor. | |
virtual | ~CdrGodunov () |
Destructor (does nothing) | |
virtual void | parseOptions () override |
Parse class options to put object in usable state. | |
virtual void | parseRuntimeOptions () override |
Parse runtime options. | |
virtual void | allocate () override |
Allocate internal data holders in class object. | |
virtual Real | computeAdvectionDt () override |
Compute the largest possible advective time step (for explicit methods) More... | |
Public Member Functions inherited from CdrMultigrid | |
CdrMultigrid () | |
Constructor. | |
virtual | ~CdrMultigrid () |
Destructor. | |
virtual void | registerOperators () override |
Register operator. | |
virtual void | preRegrid (const int a_lbase, const int a_oldFinestLevel) override |
Perform pre-regrid operations. More... | |
virtual void | computeDivJ (EBAMRCellData &a_divJ, EBAMRCellData &a_phi, const Real a_extrapDt, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final |
Compute div(J) explicitly, where J = nV - D*grad(n) More... | |
virtual void | computeDivF (EBAMRCellData &a_divF, EBAMRCellData &a_phi, const Real a_extrapDt, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final |
Compute div(v*phi) explicitly. More... | |
virtual void | computeDivD (EBAMRCellData &a_divD, EBAMRCellData &a_phi, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final |
Compute div(D*grad(phi)) explicitly. More... | |
virtual void | advanceEuler (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt) override final |
Implicit diffusion Euler advance with source term. More... | |
virtual void | advanceCrankNicholson (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt) override final |
Implicit diffusion Crank-Nicholson advance with source term. More... | |
Public Member Functions inherited from CdrSolver | |
CdrSolver () | |
Default constructor. More... | |
CdrSolver (const CdrSolver &a_other)=delete | |
Disallowed copy constructor. More... | |
CdrSolver (const CdrSolver &&a_other)=delete | |
Disallowed move constructor. More... | |
CdrSolver & | operator= (const CdrSolver &a_other)=delete |
Disallowed assignment operator. More... | |
CdrSolver & | operator= (const CdrSolver &&a_other)=delete |
Disallowed move assignement operator. More... | |
virtual | ~CdrSolver () |
Constructor. | |
void | setDefaultDomainBC () |
This sets default boundary conditions (wall type). More... | |
void | setDomainBcType (const CdrDomainBC::DomainSide a_domainSide, const CdrDomainBC::BcType a_bcType) |
Set domain bc type on domain side. More... | |
void | setDomainBcFunction (const CdrDomainBC::DomainSide a_domainSide, const CdrDomainBC::FluxFunction a_fluxFunction) |
Set domain bc function on particular domain side. More... | |
virtual void | advanceEuler (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const Real a_dt) |
Implicit diffusion Euler advance without source term. More... | |
virtual void | advanceCrankNicholson (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const Real a_dt) |
Implicit diffusion Crank-Nicholson advance without source term. More... | |
virtual void | computeDivG (EBAMRCellData &a_divG, EBAMRFluxData &a_G, const EBAMRIVData &a_ebFlux, const bool a_conservativeOnly) |
Compute div(G) where G is a general face-centered flux on face centers and EB centers. This can involve mass redistribution. More... | |
virtual void | gwnDiffusionSource (EBAMRCellData &a_noiseSource, const EBAMRCellData &a_cellPhi) |
Compute a random gaussian white noise source term. More... | |
virtual void | averageVelocityToFaces () |
Average velocities to faces. More... | |
virtual void | deallocate () |
Deallocate internal storage. More... | |
virtual void | setRealm (const std::string a_realm) |
Set the realm for this solver. More... | |
virtual RefCountedPtr< CdrSpecies > & | getSpecies () noexcept |
Get the CDR species. | |
virtual const RefCountedPtr< CdrSpecies > & | getSpecies () const noexcept |
Get the CDR species. | |
virtual void | setSpecies (const RefCountedPtr< CdrSpecies > &a_species) |
Set species. More... | |
virtual void | setComputationalGeometry (const RefCountedPtr< ComputationalGeometry > &a_computationalGeometry) |
Set computational geometry. More... | |
virtual void | setAmr (const RefCountedPtr< AmrMesh > &a_amr) |
Set the amr object. More... | |
virtual void | setPhase (phase::which_phase a_phase) |
Set phase. More... | |
virtual void | setVerbosity (const int a_verbosity) |
Set verbosity. More... | |
virtual void | setTime (const int a_step, const Real a_time, const Real a_dt) |
Set the time for this solver. More... | |
virtual void | setVelocity (const EBAMRCellData &a_velocity) |
Set velocity from data holder. More... | |
virtual void | setVelocity (const RealVect a_velocity) |
Set constant velocity. More... | |
virtual void | setVelocity (const std::function< RealVect(const RealVect a_pos)> &a_velocity) |
Set velocity using a polymorphic function. More... | |
virtual void | setDiffusionCoefficient (const EBAMRFluxData &a_diffusionCoefficient, const EBAMRIVData &a_ebDiffusionCoefficient) |
Data-based version of setting diffusion coefficients (which are stored on faces) More... | |
virtual void | setDiffusionCoefficient (const Real a_diffusionCoefficient) |
Set diffusion coefficient to be constant everywhere (they are stored on faces) More... | |
virtual void | setDiffusionCoefficient (const std::function< Real(const RealVect a_position)> &a_diffusionCoefficient) |
Polymorphic way of setting diffusion coefficients every. More... | |
virtual void | setSource (const EBAMRCellData &a_source) |
Data based version of setting source terms. More... | |
virtual void | setSource (const Real a_source) |
Set constant source terms everywhere. More... | |
virtual void | setSource (const std::function< Real(const RealVect a_position)> a_source) |
Polymorphic way of setting source terms. More... | |
virtual void | setEbFlux (const EBAMRIVData &a_ebFlux) |
Data-based version of setting the EB flux. More... | |
virtual void | setEbFlux (const Real a_ebFlux) |
Set the eb flux to a constant. More... | |
virtual void | initialData () |
Fill m_phi state with initial data from m_species. More... | |
virtual void | writePlotFile () |
Write plot file. More... | |
virtual void | writePlotData (LevelData< EBCellFAB > &a_output, int &a_icomp, const std::string a_outputRealm, const int a_level) const noexcept |
Write output data to a_output. More... | |
virtual void | regrid (const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) |
Write checkpoint data into HDF5 file. @paramo[out] a_handle HDF5 file. More... | |
virtual std::string | getRealm () const |
Get the realm where this solver is registered. More... | |
virtual std::string | getName () const |
Get solver name. More... | |
virtual Vector< std::string > | getPlotVariableNames () const |
Get output plot names. More... | |
virtual int | getNumberOfPlotVariables () const |
Get number of output fields. More... | |
virtual Real | computeDiffusionDt () |
Compute the largest possible diffusive time step (for explicit methods) More... | |
virtual Real | computeAdvectionDiffusionDt () |
Compute the largest possible diffusive time step (for explicit methods) More... | |
virtual Real | computeSourceDt (const Real a_max, const Real a_tolerance) |
Compute the largest possible source time step (for explicit methods. More... | |
virtual void | weightedUpwind (EBAMRCellData &a_weightedUpwindPhi, const int a_pow) |
Compute an upwind-weighted version of phi. More... | |
virtual Real | computeMass () |
Compute the "physical mass" in m_phi. More... | |
virtual Real | computeMass (const EBAMRCellData &a_phi, const bool a_kappaScale=true) |
Compute the "physical mass" in the input argument. More... | |
virtual Real | computeCharge () |
Compute the total charge in m_phi. More... | |
virtual bool | isDiffusive () |
Return true if the solver is diffusive and false otherwise. | |
virtual bool | isMobile () |
Return true if the solver is mobile and false otherwise. | |
virtual EBAMRCellData & | getPhi () |
Get the cell-centered phi. More... | |
virtual EBAMRCellData & | getSource () |
Get the source term. More... | |
virtual EBAMRCellData & | getCellCenteredVelocity () |
Get the cell-centered velocity. More... | |
virtual EBAMRFluxData & | getFaceCenteredVelocity () |
Get the face-centered velocities. More... | |
virtual EBAMRIVData & | getEbCenteredVelocity () |
Get the eb-centered velocities. More... | |
virtual EBAMRCellData & | getCellCenteredDiffusionCoefficient () |
Get the cell-centered diffusion coefficient. | |
virtual EBAMRFluxData & | getFaceCenteredDiffusionCoefficient () |
Get the face-centered diffusion coefficient. More... | |
virtual EBAMRIVData & | getEbCenteredDiffusionCoefficient () |
Get the EB-centered diffusion coefficient. More... | |
virtual EBAMRIVData & | getEbFlux () |
Get the eb flux data holder. More... | |
virtual EBAMRIFData & | getDomainFlux () |
Get the domain flux data holder. More... | |
virtual void | extrapolateAdvectiveFluxToEB () noexcept |
Extrapolate advective flux to EB. More... | |
virtual void | extrapolateAdvectiveFluxToEB (EBAMRIVData &a_ebFlux) const noexcept |
Extrapolate advective flux to EB. More... | |
virtual void | redistribute (EBAMRCellData &a_phi, const EBAMRIVData &a_delta) const noexcept |
Add data through redistribution into cell-centered holders. More... | |
Protected Member Functions | |
virtual void | parseSlopeLimiter () |
Parses slope limiter options. | |
virtual void | parseExtrapolateSourceTerm () |
virtual void | advectToFaces (EBAMRFluxData &a_facePhi, const EBAMRCellData &a_phi, const Real a_extrapDt) override |
Godunov face extrapolation method for advection. More... | |
Protected Member Functions inherited from CdrMultigrid | |
virtual void | setupDiffusionSolver () |
Set up diffusion solver. | |
virtual void | setupHelmholtzFactory () |
Setup the operator factory. | |
virtual void | setupMultigrid () |
Setup multigrid. | |
virtual void | setMultigridSolverCoefficients () |
Set the multigrid solver coefficients. More... | |
virtual void | resetAlphaAndBeta (const Real a_alpha, const Real a_beta) |
Reset alpha and beta-coefficients in the multigrid solvers. More... | |
virtual void | computeKappaLphi (EBAMRCellData &a_kappaLphi, const EBAMRCellData &a_phi) |
Compute kappa * L(phi) using the multigrid operator. This is different from computeDivD. More... | |
virtual void | parseMultigridSettings () |
Parse solver settings for geometric multigrid. | |
Protected Member Functions inherited from CdrSolver | |
virtual void | averageVelocityToFaces (EBAMRFluxData &a_faceVelocity, const EBAMRCellData &a_cellVelocity) |
Average cell-centered velocities to face centers. More... | |
virtual void | computeAdvectionFlux (EBAMRFluxData &a_flux, const EBAMRFluxData &a_facePhi, const EBAMRFluxData &a_faceVelocity, const bool a_addDomainFlux=true) |
Set up face-centered advection flux. More... | |
virtual void | computeAdvectionFlux (LevelData< EBFluxFAB > &a_flux, const LevelData< EBFluxFAB > &a_facePhi, const LevelData< EBFluxFAB > &a_faceVelocity, const int a_lvl) |
Set up face-centered advection flux on a grid level. More... | |
virtual void | computeDiffusionFlux (EBAMRFluxData &a_flux, const EBAMRCellData &a_phi, const bool a_addDomainFlux) |
Compute the face-centered diffusion flux. More... | |
virtual void | computeDiffusionFlux (LevelData< EBFluxFAB > &a_flux, const LevelData< EBCellFAB > &a_phi, const int a_lvl) |
Compute the face-centered diffusion flux. More... | |
virtual void | computeAdvectionDiffusionFlux (EBAMRFluxData &a_flux, const EBAMRCellData &a_cellStates, const EBAMRFluxData &a_faceStates, const EBAMRFluxData &a_faceVelocities, const EBAMRFluxData &a_faceDiffCo, const bool a_addDomainFlux) |
Compute the full advection-diffusion flux. This assumes that the solver is mobile and diffusive. More... | |
virtual void | resetDomainFlux (EBAMRFluxData &a_flux) |
Set flux to zero on domain boundaries. More... | |
virtual void | fillDomainFlux (EBAMRFluxData &a_flux) |
Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux. More... | |
virtual void | fillDomainFlux (LevelData< EBFluxFAB > &a_flux, const int a_level) |
Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux. More... | |
virtual void | conservativeDivergenceNoKappaDivision (EBAMRCellData &a_conservativeDivergence, EBAMRFluxData &a_flux, const EBAMRIVData &a_ebFlux) |
Compute conservative divergence from fluxes. More... | |
virtual void | nonConservativeDivergence (EBAMRIVData &a_nonConservativeDivergence, const EBAMRCellData &a_divG) |
Compute the non-conservative divergence. More... | |
virtual void | hybridDivergence (EBAMRCellData &a_hybridDivergence, EBAMRIVData &a_massDifference, const EBAMRIVData &a_nonConservativeDivergence) |
Use the non-conservative divergence to make the conservative divergence hold the hybrid divergence. More... | |
virtual void | hybridDivergence (LevelData< EBCellFAB > &a_hybridDivergence, LevelData< BaseIVFAB< Real >> &a_massDifference, const LevelData< BaseIVFAB< Real >> &a_nonConservativeDivergence, const int a_lvl) |
Make the hybrid divergence. On the way in, a_hybridDivergence must hold the conservative divergence. More... | |
virtual void | conservativeDivergenceRegular (LevelData< EBCellFAB > &a_divJ, const LevelData< EBFluxFAB > &a_flux, const int a_lvl) |
Compute the conservative divergence over regular cells. More... | |
virtual void | interpolateFluxToFaceCentroids (EBAMRFluxData &a_flux) |
Interpolate flux to centroids. More... | |
virtual void | interpolateFluxToFaceCentroids (LevelData< EBFluxFAB > &a_flux, const int a_lvl) |
Interpolate flux to centroids. More... | |
virtual void | computeDivergenceIrregular (LevelData< EBCellFAB > &a_divG, const LevelData< EBFluxFAB > &a_centroidFluxes, const LevelData< BaseIVFAB< Real >> &a_ebFlux, const int a_lvl) |
Compute conservative divergence on irregular cells (not kappa divided) More... | |
virtual void | initialDataDistribution () |
Fill initial data from a distribution function. More... | |
virtual void | initialDataParticles () |
Fill initial data from particles. More... | |
virtual void | defineInterpolationStencils () |
Define stencils for doing face-centered to face-centroid-centered states. More... | |
virtual void | parseDomainBc () |
Parses domain BC options. | |
virtual void | parsePlotVariables () |
Parses plot variables. | |
virtual void | parsePlotMode () |
Parse plot mode. | |
virtual void | parseDivergenceComputation () |
Parse the conservation. | |
virtual void | parseRegridSlopes () |
Parse slope regrid. | |
virtual std::string | makeBcString (const int a_dir, const Side::LoHiSide a_side) const |
Shortcut for making a boundary condition string. More... | |
virtual void | fillGwn (EBAMRFluxData &a_noise, const Real a_sigma) |
Gaussian noise field. More... | |
virtual void | smoothHeavisideFaces (EBAMRFluxData &a_facePhi, const EBAMRCellData &a_cellPhi) |
Use Heaviside smoothing for computing face-centered states. More... | |
virtual void | writeData (LevelData< EBCellFAB > &a_output, int &a_comp, const EBAMRCellData &a_data, const std::string a_outputRealm, const int a_level, const bool a_interpToCentroids, const bool a_interpGhost) const noexcept |
Write data to output. Convenience function. More... | |
Protected Attributes | |
Vector< RefCountedPtr< EBAdvectLevelIntegrator > > | m_levelAdvect |
Advection object – this implements the Graves/Trebotich discretization. | |
bool | m_limitSlopes |
If true, slopes are limited (always use limiting). | |
bool | m_extrapolateSourceTerm |
Turn on/off source terms when time-extrapolating. | |
Protected Attributes inherited from CdrMultigrid | |
EBHelmholtzOp::Smoother | m_smoother |
Relaxation type for gmg. | |
MultigridType | m_multigridType |
GMG multigrid type. | |
RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > | m_multigridSolver |
Geometric multigrid solver. | |
RefCountedPtr< EBHelmholtzOpFactory > | m_helmholtzOpFactory |
Operator factory. | |
BiCGStabSolver< LevelData< EBCellFAB > > | m_bicgstab |
Conjugate gradient solver bottom MG level. | |
EBSimpleSolver | m_simpleSolver |
Simple solver. | |
GMRESSolver< LevelData< EBCellFAB > > | m_gmres |
GMRES solver. | |
bool | m_hasMultigridSolver |
Is multigrid solver defined or not? | |
EBAMRCellData | m_helmAcoef |
Storage for Helmholtz a-coefficient. Always 1. | |
EBAMRCellData | m_residual |
Storage for multigrid residual. | |
int | m_multigridVerbosity |
Verbosity for geometric multigrid. | |
int | m_multigridPreSmooth |
Number of smoothings before averaging. | |
int | m_multigridPostSmooth |
Number of smoothings before averaging. | |
int | m_multigridBottomSmooth |
Number of smoothing before bottom solver. | |
int | m_multigridMaxIterations |
Maximum number of iterations. | |
int | m_multigridMinIterations |
Minimum number of iterations. | |
BottomSolverType | m_bottomSolverType |
Bottom solver type. | |
int | m_numSmoothingsForSimpleSolver |
Number of smoothing for bottom solver. | |
int | m_minCellsBottom |
Set bottom drop depth. | |
Real | m_multigridExitTolerance |
Multigrid tolerance. | |
Real | m_multigridExitHang |
Multigrid hand. | |
Protected Attributes inherited from CdrSolver | |
RefCountedPtr< CdrSpecies > | m_species |
Species through which e.g. mobility/diffusion and initial conditions is passed. | |
RefCountedPtr< ComputationalGeometry > | m_computationalGeometry |
Computational geometry. | |
RefCountedPtr< AmrMesh > | m_amr |
AMR; needed for grid stuff. | |
Vector< RefCountedPtr< LayoutData< BaseIFFAB< FaceStencil > > > > | m_interpStencils [SpaceDim] |
Stencils for interpolating face-centered fluxes to face centroids. | |
phase::which_phase | m_phase |
Phase. | |
std::string | m_name |
Solver name. | |
std::string | m_className |
Class name. More... | |
std::string | m_realm |
Realm where this solver is registered. | |
CdrDomainBC | m_domainBC |
Domain BCs. | |
EBAMRCellData | m_phi |
Cell-centered data (i.e. the advected-diffused quantity) | |
EBAMRCellData | m_source |
Source term. | |
EBAMRCellData | m_cellVelocity |
Cell-centered velocities. More... | |
EBAMRFluxData | m_faceStates |
Holder for face centered states. More... | |
EBAMRFluxData | m_faceVelocity |
Face-centered velocities (only normal components) More... | |
EBAMRIVData | m_nonConservativeDivG |
Scratch storage for the non-conservative divergence. | |
EBAMRIVData | m_massDifference |
Scratch storage for the mass difference. | |
EBAMRIVData | m_ebZero |
Scratch storage for the EB flux. | |
EBAMRCellData | m_cachePhi |
Cached state vector for regrid. | |
EBAMRCellData | m_cacheSource |
Cached source term for regrids. | |
EBAMRIVData | m_ebVelocity |
EB-centered velocities. | |
EBAMRIVData | m_ebFlux |
Flux through the embedded boundary. | |
EBAMRIFData | m_domainFlux |
Domain flux. | |
EBAMRCellData | m_cellCenteredDiffusionCoefficient |
Diffusion coefficients on cell centers. More... | |
EBAMRFluxData | m_faceCenteredDiffusionCoefficient |
Diffusion coefficients on face centers. More... | |
EBAMRIVData | m_ebCenteredDiffusionCoefficient |
Diffusion coefficients on EB faces. More... | |
std::map< CdrDomainBC::DomainSide, CdrDomainBC::FluxFunction > | m_domainFluxFunctions |
Domain flux functions. More... | |
int | m_verbosity |
Solver verbosity. More... | |
int | m_timeStep |
Time step. | |
Real | m_time |
Current time. | |
Real | m_dt |
Last time step increment. | |
Redistribution | m_whichRedistribution |
Which type of redistribution to use. | |
bool | m_blendConservation |
Flag for blending the hybrid divergence. | |
bool | m_isDiffusive |
Is the solver diffusive or not. | |
bool | m_isMobile |
Solve for advection/convection or not. | |
bool | m_plotPhi |
If true, m_phi is added to plot files. | |
bool | m_plotVelocity |
Output velocities. | |
bool | m_plotDiffusionCoefficient |
Output diffusion coefficients. | |
bool | m_plotEbFlux |
Output EB fluxes. | |
bool | m_plotSource |
Output source term. | |
bool | m_plotNumbers |
Plot numbers or densities. | |
bool | m_regridSlopes |
Use slopes when regridding. | |
int | m_seed |
RNG seed. | |
Additional Inherited Members | |
Protected Types inherited from CdrMultigrid | |
enum class | BottomSolverType { Simple , BiCGStab , GMRES } |
Enum class for supported bottom solvers in multigrid. | |
enum class | MultigridType { VCycle , WCycle } |
Enum for multigrid cycle types. | |
Protected Types inherited from CdrSolver | |
enum class | Redistribution { VolumeWeighted , None } |
Redistribution method. | |
Static Protected Attributes inherited from CdrSolver | |
static constexpr int | m_comp = 0 |
Component number in data holder. | |
static constexpr int | m_nComp = 1 |
Number of components that this solver solves for. | |
Godunov implementation for advection.
The CdrGodunov class perform Godunov-type advection with time-extrapolated slope-limiting to faces. Note that this class implements CdrMultigrid which is the base class that implements second order diffusion (either explicit or implicit). The rationale for doing time-extrapolation ala Godunov is that one can use centered differencing phi^(k+1) = phi^k - dt*div[v^(k+1/2) * phi^(k+1/2)] for a second order accurate discretization in time. This class uses Chombo functionality from EBAdvectLevelIntegrator which uses the Graves/Trebotich discretization in "An adaptive finite volume method for the incompressible Navier–Stokes equations in complex geometries" (DOI: 10.2140/camcos.2015.10.43). Note that if source terms are involved the extrapolation is more involved.
|
overrideprotectedvirtual |
Godunov face extrapolation method for advection.
[out] | a_facePhi | Phi on face centers |
[in] | a_cellPhi | Phi on cell centers |
[in] | a_extrapDt | Time centering (i.e. extrapolation) of the face-centered states. |
This functions calls Chombo code which implements the Graves/Trebotich discretization.
Implements CdrMultigrid.
|
overridevirtual |
Compute the largest possible advective time step (for explicit methods)
This computes dt = dx/max(|vx|,|vy|,|vz|), minimized over all grid levels and patches.
Reimplemented from CdrSolver.
|
protectedvirtual |
®brief Parse whether or not to include a source-term in the Godunov extrapolation.