chombo-discharge
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
CdrSolver Class Referenceabstract

Base class for solving convection-diffusion-reaction equations. More...

#include <CD_CdrSolver.H>

Inheritance diagram for CdrSolver:
Inheritance graph
[legend]
Collaboration diagram for CdrSolver:
Collaboration graph
[legend]

Public Member Functions

 CdrSolver ()
 Default constructor. More...
 
 CdrSolver (const CdrSolver &a_other)=delete
 Disallowed copy constructor. More...
 
 CdrSolver (const CdrSolver &&a_other)=delete
 Disallowed move constructor. More...
 
CdrSolveroperator= (const CdrSolver &a_other)=delete
 Disallowed assignment operator. More...
 
CdrSolveroperator= (const CdrSolver &&a_other)=delete
 Disallowed move assignement operator. More...
 
virtual ~CdrSolver ()
 Constructor.
 
virtual void parseOptions ()=0
 Parse class options.
 
virtual void parseRuntimeOptions ()=0
 Parse runtime options.
 
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 advanceEuler (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt)=0
 Implicit diffusion Euler advance with 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 advanceCrankNicholson (EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt)=0
 Implicit diffusion Crank-Nicholson advance with source term. 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)=0
 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)=0
 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)=0
 Compute div(D*grad(phi)) explicitly. 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 preRegrid (const int a_lbase, const int a_oldFinestLevel)
 Perform pre-regrid operations. More...
 
virtual void deallocate ()
 Deallocate internal storage. More...
 
virtual void registerOperators ()
 Register operators for AMR operations. 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 allocate ()
 Allocate internal storage. 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 computeAdvectionDt ()
 Compute the largest possible diffusive time step (for explicit methods) 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 EBAMRCellDatagetPhi ()
 Get the cell-centered phi. More...
 
virtual EBAMRCellDatagetSource ()
 Get the source term. More...
 
virtual EBAMRCellDatagetCellCenteredVelocity ()
 Get the cell-centered velocity. More...
 
virtual EBAMRFluxDatagetFaceCenteredVelocity ()
 Get the face-centered velocities. More...
 
virtual EBAMRIVDatagetEbCenteredVelocity ()
 Get the eb-centered velocities. More...
 
virtual EBAMRCellDatagetCellCenteredDiffusionCoefficient ()
 Get the cell-centered diffusion coefficient.
 
virtual EBAMRFluxDatagetFaceCenteredDiffusionCoefficient ()
 Get the face-centered diffusion coefficient. More...
 
virtual EBAMRIVDatagetEbCenteredDiffusionCoefficient ()
 Get the EB-centered diffusion coefficient. More...
 
virtual EBAMRIVDatagetEbFlux ()
 Get the eb flux data holder. More...
 
virtual EBAMRIFDatagetDomainFlux ()
 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 Types

enum class  Redistribution { VolumeWeighted , None }
 Redistribution method.
 

Protected Member Functions

virtual void averageVelocityToFaces (EBAMRFluxData &a_faceVelocity, const EBAMRCellData &a_cellVelocity)
 Average cell-centered velocities to face centers. More...
 
virtual void advectToFaces (EBAMRFluxData &a_facePhi, const EBAMRCellData &a_phi, const Real a_extrapDt)=0
 Advection-only extrapolation to faces. 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

RefCountedPtr< CdrSpeciesm_species
 Species through which e.g. mobility/diffusion and initial conditions is passed.
 
RefCountedPtr< ComputationalGeometrym_computationalGeometry
 Computational geometry.
 
RefCountedPtr< AmrMeshm_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::FluxFunctionm_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.
 

Static Protected Attributes

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.
 

Detailed Description

Base class for solving convection-diffusion-reaction equations.

On embedded boundaries, users must always set the flux directly (but possibly via a function, if desired).

On domain edges/faces, the user can select between various ways of setting the domain boundary conditions (e.g. wall, data-based, function-based, outflow).

Constructor & Destructor Documentation

◆ CdrSolver() [1/3]

CdrSolver::CdrSolver ( )

Default constructor.

Note
This sets the realm to primal and populates the domain flux functions (with zero-fluxes, i.e. wall BCs)

◆ CdrSolver() [2/3]

CdrSolver::CdrSolver ( const CdrSolver a_other)
delete

Disallowed copy constructor.

Parameters
[in]a_otherOther solver

◆ CdrSolver() [3/3]

CdrSolver::CdrSolver ( const CdrSolver &&  a_other)
delete

Disallowed move constructor.

Parameters
[in]a_otherOther solver

Member Function Documentation

◆ advanceCrankNicholson() [1/2]

virtual void CdrSolver::advanceCrankNicholson ( EBAMRCellData a_newPhi,
const EBAMRCellData a_oldPhi,
const EBAMRCellData a_source,
const Real  a_dt 
)
pure virtual

Implicit diffusion Crank-Nicholson advance with source term.

Parameters
[in,out]a_newPhiSolution at time t + dt
[in]a_oldPhiSolution at time t
[in]a_sourceSource term.
[in]a_dtTime step

Implemented in CdrMultigrid.

◆ advanceCrankNicholson() [2/2]

void CdrSolver::advanceCrankNicholson ( EBAMRCellData a_newPhi,
const EBAMRCellData a_oldPhi,
const Real  a_dt 
)
virtual

Implicit diffusion Crank-Nicholson advance without source term.

Parameters
[in,out]a_newPhiSolution at time t + dt
[in]a_oldPhiSolution at time t
[in]a_dtTime step
Note
This calls the other version with a zero source term.

◆ advanceEuler() [1/2]

virtual void CdrSolver::advanceEuler ( EBAMRCellData a_newPhi,
const EBAMRCellData a_oldPhi,
const EBAMRCellData a_source,
const Real  a_dt 
)
pure virtual

Implicit diffusion Euler advance with source term.

Parameters
[in,out]a_newPhiSolution at time t + dt
[in]a_oldPhiSolution at time t
[in]a_sourceSource term.
[in]a_dtTime step
Note
For purely implicit Euler the source term should be centered at t+dt (otherwise it's an implicit-explicit method)

This solves the implicit diffusion equation equation a_newPhi - a_oldPhi = dt*Laplacian(a_newPhi) + dt*a_source.

Implemented in CdrMultigrid.

◆ advanceEuler() [2/2]

void CdrSolver::advanceEuler ( EBAMRCellData a_newPhi,
const EBAMRCellData a_oldPhi,
const Real  a_dt 
)
virtual

Implicit diffusion Euler advance without source term.

Parameters
[in,out]a_newPhiSolution at time t + dt
[in]a_oldPhiSolution at time t
[in]a_dtTime step
Note
This calls the other version with a zero source term.

◆ advectToFaces()

virtual void CdrSolver::advectToFaces ( EBAMRFluxData a_facePhi,
const EBAMRCellData a_phi,
const Real  a_extrapDt 
)
protectedpure virtual

Advection-only extrapolation to faces.

Parameters
[out]a_facePhiPhi on faces
[in]a_phiPhi on cell center
[in]a_extrapDtTime centering/extrapolation (if the advective integrator can do it)

Implemented in CdrMultigrid, CdrGodunov, and CdrCTU.

◆ allocate()

void CdrSolver::allocate ( )
virtual

Allocate internal storage.

Note
This allocates a bunch of storage – if the solver is not diffusive or mobile then we do our best to trim memory.

Reimplemented in CdrMultigrid, and CdrGodunov.

◆ averageVelocityToFaces() [1/2]

void CdrSolver::averageVelocityToFaces ( )
virtual

Average velocities to faces.

Note
This takes the cell-centered velocities in m_cellVelocity and averages them onto m_faceVelocity. Only the velocity component normal to the face is stored in m_faceVelocity.

◆ averageVelocityToFaces() [2/2]

void CdrSolver::averageVelocityToFaces ( EBAMRFluxData a_faceVelocity,
const EBAMRCellData a_cellVelocity 
)
protectedvirtual

Average cell-centered velocities to face centers.

Parameters
[out]a_faceVelocityFace-centered velocities
[in]a_cellVelocityCell-centered velocities
Note
This sets the velocity component normal to the face to be the arithmetic average of the cell center velocities on the low/high side of the face.
Virtual because some solvers might want to do this differently.

◆ computeAdvectionDiffusionDt()

Real CdrSolver::computeAdvectionDiffusionDt ( )
virtual

Compute the largest possible diffusive time step (for explicit methods)

This computes dt = 1/[ dx/(|vx|+|vy|+|vz|) + (dx*dx)/(2*D*d) ] and minimizes the result over all grid levels and patches.

Note
This is the appropriate time step routine for explicit advection-diffusion solvers.

◆ computeAdvectionDiffusionFlux()

void CdrSolver::computeAdvectionDiffusionFlux ( EBAMRFluxData a_flux,
const EBAMRCellData a_cellStates,
const EBAMRFluxData a_faceStates,
const EBAMRFluxData a_faceVelocities,
const EBAMRFluxData a_faceDiffCo,
const bool  a_addDomainFlux 
)
protectedvirtual

Compute the full advection-diffusion flux. This assumes that the solver is mobile and diffusive.

Parameters
[out]a_fluxFace-centered flux.
[in]a_cellStatesCell-centered states.
[in]a_faceVelocitiesFace-centered velocities.
[in]a_faceStatesFace-centered (advected) states.
[in]a_faceDiffCoFace-centered diffusion coefficients.
[in]a_addDomainFluxAdd domain flux or not.

◆ computeAdvectionDt()

Real CdrSolver::computeAdvectionDt ( )
virtual

Compute the largest possible diffusive time step (for explicit methods)

This computes dt = dx/sum(|vx| + |vy| + |vz|), minimized over all grid levels and patches.

Note
This is the appropriate time step routine for explicit advection solvers.

Reimplemented in CdrGodunov, and CdrCTU.

◆ computeAdvectionFlux() [1/2]

void CdrSolver::computeAdvectionFlux ( EBAMRFluxData a_flux,
const EBAMRFluxData a_facePhi,
const EBAMRFluxData a_faceVelocity,
const bool  a_addDomainFlux = true 
)
protectedvirtual

Set up face-centered advection flux.

Parameters
[out]a_fluxFace-centered fluxes
[in]a_facePhiFace-centered state
[in]a_faceVelocityFace-centered velocities
[in]a_addDomainFluxAdd the domain flux
Note
This computes flux = phi*vel on all face centers. If a_addDomainFlux is true, we enforce BC fluxes on the domain sides. Otherwise, the flux on the domain sides is set to zero.

◆ computeAdvectionFlux() [2/2]

void CdrSolver::computeAdvectionFlux ( LevelData< EBFluxFAB > &  a_flux,
const LevelData< EBFluxFAB > &  a_facePhi,
const LevelData< EBFluxFAB > &  a_faceVelocity,
const int  a_lvl 
)
protectedvirtual

Set up face-centered advection flux on a grid level.

Parameters
[out]a_fluxFace-centered fluxes
[in]a_facePhiFace-centered state
[in]a_faceVelocityFace-centered velocities
[in]a_domainFluxDomain flux
[in]a_levelGrid level
Note
This computes flux = v*phi on all face centers. If a_addDomainFlux is true, we enforce BC fluxes on the domain sides. Otherwise, the flux on the domain sides is set to zero.

◆ computeCharge()

Real CdrSolver::computeCharge ( )
virtual

Compute the total charge in m_phi.

Calls computeMass and multiplies by charge

◆ computeDiffusionDt()

Real CdrSolver::computeDiffusionDt ( )
virtual

Compute the largest possible diffusive time step (for explicit methods)

This computes dt = (dx*dx)/(2*D*d) where D is the diffusion coefficient. The result is minimized over all grid levels and patches.

Note
This is the appropriate time step routine for explicit diffusion solvers.

◆ computeDiffusionFlux() [1/2]

void CdrSolver::computeDiffusionFlux ( EBAMRFluxData a_flux,
const EBAMRCellData a_phi,
const bool  a_addDomainFlux 
)
protectedvirtual

Compute the face-centered diffusion flux.

Parameters
[out]a_fluxFace-centered flux
[in]a_phiCell centered phi
[in]a_addDomainFluxAdd domain flux or not
Note
This computes flux = D*grad(phi) on all face centers. If a_addDomainFlux is true, we enforce BC fluxes on the domain sides. Otherwise, the flux on the domain sides is set to zero.

◆ computeDiffusionFlux() [2/2]

void CdrSolver::computeDiffusionFlux ( LevelData< EBFluxFAB > &  a_flux,
const LevelData< EBCellFAB > &  a_phi,
const int  a_lvl 
)
protectedvirtual

Compute the face-centered diffusion flux.

Parameters
[out]a_fluxFace-centered flux
[in]a_phiCell centered phi
[in]a_lvlGrid level
Note
This computes flux = D*grad(phi) on all face centers. The flux is set to zero on domain edges/faces

◆ computeDivD()

virtual void CdrSolver::computeDivD ( EBAMRCellData a_divD,
EBAMRCellData a_phi,
const bool  a_conservativeOnly,
const bool  a_ebFlux,
const bool  a_domainFlux 
)
pure virtual

Compute div(D*grad(phi)) explicitly.

Parameters
[out]a_divFDivergence term, i.e. finite volume approximation to Div(D*Grad(phi)).
[in]a_phiCell-centered state
[in]a_domainBcFlag for setting domain boundary conditions
[in]a_conservativeOnlyIf true, we compute div(D) = 1/dx*sum(fluxes), which does not involve redistribution.
[in]a_useEbFluxIf true, the embedded boundary flux will be injected and included in div(D)
[in]a_domainFluxIf true, the domain flux will injected and included in div(D)
Note
a_phi is non-const because ghost cells will be interpolated in this routine. Valid data in a_phi is not touched.

Implemented in CdrMultigrid.

◆ computeDivergenceIrregular()

void CdrSolver::computeDivergenceIrregular ( LevelData< EBCellFAB > &  a_divG,
const LevelData< EBFluxFAB > &  a_centroidFluxes,
const LevelData< BaseIVFAB< Real >> &  a_ebFlux,
const int  a_lvl 
)
protectedvirtual

Compute conservative divergence on irregular cells (not kappa divided)

Parameters
[out]a_divGConservative divergence (not kappa divided)
[out]a_centroidFluxFace centroid centered fluxes
[in]a_ebFluxEB flux
[in]a_lvlGrid level
Note
This overwrites the result of divG in cut-cells, using the eb flux and face centroid fluxes to compute kappa*div(G) = sum(fluxes)/dx

◆ computeDivF()

virtual void CdrSolver::computeDivF ( EBAMRCellData a_divF,
EBAMRCellData a_phi,
const Real  a_extrapDt,
const bool  a_conservativeOnly,
const bool  a_ebFlux,
const bool  a_domainFlux 
)
pure virtual

Compute div(v*phi) explicitly.

Parameters
[out]a_divFDivergence term, i.e. finite volume approximation to Div(v*phi), including redistribution magic.
[in]a_phiCell-centered state
[in]a_extrapDtExtrapolation in time, i.e. shifting of div(F) towards e.g. half time step. Only affects the advective term.
[in]a_conservativeOnlyIf true, we compute div(F)=1/dx*sum(fluxes), which does not involve redistribution.
[in]a_domainBcHow to set domain fluxes
[in]a_ebFluxIf true, the embedded boundary flux will injected be included in div(F)
[in]a_domainFluxIf true, the domain flux will injected and included in div(F)
Note
a_phi is non-const because ghost cells will be interpolated in this routine. Valid data in a_phi is not touched.

Implemented in CdrMultigrid.

◆ computeDivG()

void CdrSolver::computeDivG ( EBAMRCellData a_divG,
EBAMRFluxData a_G,
const EBAMRIVData a_ebFlux,
const bool  a_conservativeOnly 
)
virtual

Compute div(G) where G is a general face-centered flux on face centers and EB centers. This can involve mass redistribution.

Parameters
[in]a_divGdiv(G) or kappa*div(G).
[in,out]a_GVector field which contains face-centered fluxes on input. Contains face-centroid fluxes on output.
[in]a_ebFluxFlux on the EB centroids
[in]a_conservativeOnlyIf true, we compute div(G)=1/dx*sum(fluxes), which does not involve redistribution.

◆ computeDivJ()

virtual void CdrSolver::computeDivJ ( EBAMRCellData a_divJ,
EBAMRCellData a_phi,
const Real  a_extrapDt,
const bool  a_conservativeOnly,
const bool  a_ebFlux,
const bool  a_domainFlux 
)
pure virtual

Compute div(J) explicitly, where J = nV - D*grad(n)

Parameters
[out]a_divJDivergence term, i.e. finite volume approximation to
[in]a_phiCell-centered state
[in]a_extrapDtExtrapolation in time, i.e. shifting of div(J) towards e.g. half time step. Only affects the advective term.
[in]a_conservativeOnlyIf true, we compute div(J) = 1/dx*sum(fluxes), which does not involve redistribution.
[in]a_ebFluxIf true, the embedded boundary flux will injected and included in div(J)
[in]a_domainFluxIf true, the domain flux will injected and included in div(J)
Note
a_phi is non-const because ghost cells will be re-filled

Implemented in CdrMultigrid.

◆ computeMass() [1/2]

Real CdrSolver::computeMass ( )
virtual

Compute the "physical mass" in m_phi.

This conservatively coarsens the solution and computes on the coarsest level only

◆ computeMass() [2/2]

Real CdrSolver::computeMass ( const EBAMRCellData a_phi,
const bool  a_kappaScale = true 
)
virtual

Compute the "physical mass" in the input argument.

This conservatively coarsens the solution and computes on the coarsest level only

Parameters
[in]a_phiCell-centered data.
[in]a_kappaScaleMultiply by kappa in irregular cells.

◆ computeSourceDt()

Real CdrSolver::computeSourceDt ( const Real  a_max,
const Real  a_tolerance 
)
virtual

Compute the largest possible source time step (for explicit methods.

Parameters
[in]a_maxMaximum value of m_phi
[in]a_toleranceTolerance

This computes dt = phi/source, but only for cells where phi lies within a_tolerance*a_max (and source > 0)

Todo:
An old routine which we could (probably) remove.

◆ conservativeDivergenceNoKappaDivision()

void CdrSolver::conservativeDivergenceNoKappaDivision ( EBAMRCellData a_conservativeDivergence,
EBAMRFluxData a_flux,
const EBAMRIVData a_ebFlux 
)
protectedvirtual

Compute conservative divergence from fluxes.

Parameters
[out]a_conservativeDivergenceConservative divergence computed as div(F) using finite volumes. Not divided by kappa.
[in]a_fluxFace-centered fluxes. Includes domain fluxes.
[in]a_ebFluxEB flux.
Note
This computes the conservative divergence kappa*div(F) = sum(fluxes).

◆ conservativeDivergenceRegular()

void CdrSolver::conservativeDivergenceRegular ( LevelData< EBCellFAB > &  a_divJ,
const LevelData< EBFluxFAB > &  a_flux,
const int  a_lvl 
)
protectedvirtual

Compute the conservative divergence over regular cells.

Parameters
[out]a_divJConservative divergence. Not kappa-divided.
[in]a_fluxFace-centered fluxes
[in]a_lvlGrid level
Note
This computes the finite volume approximation to div(J) in regular cells – cut-cells (in a_divJ are set to zero)

◆ deallocate()

void CdrSolver::deallocate ( )
virtual

Deallocate internal storage.

Note
This deallocates a bunch of storage which is not needed (during regrids). It can be used to trim memory during Berger-Rigoutsous regrids (which EAT memory).

◆ defineInterpolationStencils()

void CdrSolver::defineInterpolationStencils ( )
protectedvirtual

Define stencils for doing face-centered to face-centroid-centered states.

Note
This computes standard finite-difference stencils for interpolating from face centers to face centroids.

◆ extrapolateAdvectiveFluxToEB() [1/2]

void CdrSolver::extrapolateAdvectiveFluxToEB ( )
virtualnoexcept

Extrapolate advective flux to EB.

This does not use any slope limiting (derived classes might want to redo this routine). The returned flux is equal to -n.(phi * v) where n is the normal vector pointing into the solution domain.

Note
Calls the other version with m_ebFlux as argument

◆ extrapolateAdvectiveFluxToEB() [2/2]

void CdrSolver::extrapolateAdvectiveFluxToEB ( EBAMRIVData a_ebFlux) const
virtualnoexcept

Extrapolate advective flux to EB.

This does not use any slope limiting (derived classes might want to redo this routine). The returned flux is equal to -n.(phi * v) where n is the normal vector pointing into the solution domain.

Parameters
[out]a_ebFluxProjected flux on the EB.

◆ fillDomainFlux() [1/2]

void CdrSolver::fillDomainFlux ( EBAMRFluxData a_flux)
protectedvirtual

Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux.

Parameters
[in,out]a_fluxFlux to be modified.

◆ fillDomainFlux() [2/2]

void CdrSolver::fillDomainFlux ( LevelData< EBFluxFAB > &  a_flux,
const int  a_level 
)
protectedvirtual

Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux.

Parameters
[in,out]a_fluxFlux to be modified.
[in]a_levelGrid level

◆ fillGwn()

void CdrSolver::fillGwn ( EBAMRFluxData a_noise,
const Real  a_sigma 
)
protectedvirtual

Gaussian noise field.

Parameters
[out]a_noiseGaussian white nosie
[out]a_sigmaStandard deviation

◆ getCellCenteredVelocity()

EBAMRCellData & CdrSolver::getCellCenteredVelocity ( )
virtual

Get the cell-centered velocity.

Returns
m_cellVelocity

◆ getDomainFlux()

EBAMRIFData & CdrSolver::getDomainFlux ( )
virtual

Get the domain flux data holder.

Returns
m_domainFlux

◆ getEbCenteredDiffusionCoefficient()

EBAMRIVData & CdrSolver::getEbCenteredDiffusionCoefficient ( )
virtual

Get the EB-centered diffusion coefficient.

Returns
m_ebCentereDiffusionCoefficient

◆ getEbCenteredVelocity()

EBAMRIVData & CdrSolver::getEbCenteredVelocity ( )
virtual

Get the eb-centered velocities.

Returns
m_ebVelocity

◆ getEbFlux()

EBAMRIVData & CdrSolver::getEbFlux ( )
virtual

Get the eb flux data holder.

Returns
m_ebFlux

◆ getFaceCenteredDiffusionCoefficient()

EBAMRFluxData & CdrSolver::getFaceCenteredDiffusionCoefficient ( )
virtual

Get the face-centered diffusion coefficient.

Returns
m_faceCentereDiffusionCoefficient

◆ getFaceCenteredVelocity()

EBAMRFluxData & CdrSolver::getFaceCenteredVelocity ( )
virtual

Get the face-centered velocities.

Returns
m_faceVelocity

◆ getName()

std::string CdrSolver::getName ( ) const
virtual

Get solver name.

Note
Not necessarily equal to class name (we can have many CdrSolvers instantiated).

◆ getNumberOfPlotVariables()

int CdrSolver::getNumberOfPlotVariables ( ) const
virtual

Get number of output fields.

Returns
Returns number of plot variables include in writePlotFile() and writePlotData()

◆ getPhi()

EBAMRCellData & CdrSolver::getPhi ( )
virtual

Get the cell-centered phi.

Returns
m_phi

◆ getPlotVariableNames()

Vector< std::string > CdrSolver::getPlotVariableNames ( ) const
virtual

Get output plot names.

Returns
Return a list of plot variable names.

◆ getRealm()

std::string CdrSolver::getRealm ( ) const
virtual

Get the realm where this solver is registered.

Returns
Returns realm name.

◆ getSource()

EBAMRCellData & CdrSolver::getSource ( )
virtual

Get the source term.

Returns
m_source

◆ gwnDiffusionSource()

void CdrSolver::gwnDiffusionSource ( EBAMRCellData a_noiseSource,
const EBAMRCellData a_cellPhi 
)
virtual

Compute a random gaussian white noise source term.

Parameters
[out]a_noiseSourceSource term
[in]a_cellPhiCell-centered states
Note
Used e.g. for fluctuating hydrodynamics.

◆ hybridDivergence() [1/2]

void CdrSolver::hybridDivergence ( EBAMRCellData a_hybridDivergence,
EBAMRIVData a_massDifference,
const EBAMRIVData a_nonConservativeDivergence 
)
protectedvirtual

Use the non-conservative divergence to make the conservative divergence hold the hybrid divergence.

Parameters
[in,out]a_hybridDivergenceOn input, contains the conservative divergence (without kappa division). Contains hybrid divergence on output.
[out]a_massDifferenceMass difference between updating with hybrid divergence and true divergence
[in]a_nonConservativeDivergenceNon-conservative divergence.

◆ hybridDivergence() [2/2]

void CdrSolver::hybridDivergence ( LevelData< EBCellFAB > &  a_hybridDivergence,
LevelData< BaseIVFAB< Real >> &  a_massDifference,
const LevelData< BaseIVFAB< Real >> &  a_nonConservativeDivergence,
const int  a_lvl 
)
protectedvirtual

Make the hybrid divergence. On the way in, a_hybridDivergence must hold the conservative divergence.

Parameters
[in,out]a_hybridDivergenceOn input, contains the conservative divergence (without kappa division). Contains hybrid divergence on output.
[out]a_massDifferenceMass difference between updating with hybrid divergence and true divergence
[in]a_nonConservativeDivergenceNon-conservative divergence.
[in]a_lvlGrid level.

◆ initialData()

void CdrSolver::initialData ( )
virtual

Fill m_phi state with initial data from m_species.

Note
This will call initialDataParticles() and initialDataFunction().

◆ initialDataDistribution()

void CdrSolver::initialDataDistribution ( )
protectedvirtual

Fill initial data from a distribution function.

Note
This increments m_phi on the mesh with data from m_species->initialData().

◆ initialDataParticles()

void CdrSolver::initialDataParticles ( )
protectedvirtual

Fill initial data from particles.

Note
This will set m_phi to zero before depositing the particles (which are deposited with an NGP scheme, ignoring modifactions near cut-cells)

◆ interpolateFluxToFaceCentroids() [1/2]

void CdrSolver::interpolateFluxToFaceCentroids ( EBAMRFluxData a_flux)
protectedvirtual

Interpolate flux to centroids.

Parameters
[in,out]a_fluxOn input, contains centered fluxes. Output contains centroid fluxes

This interpolates face fluxes to centroids.

◆ interpolateFluxToFaceCentroids() [2/2]

void CdrSolver::interpolateFluxToFaceCentroids ( LevelData< EBFluxFAB > &  a_flux,
const int  a_lvl 
)
protectedvirtual

Interpolate flux to centroids.

Parameters
[in,out]a_fluxOn input, contains centered fluxes. Output contains centroid fluxes
[in]a_lvlGrid level

This interpolates face fluxes to centroids.

◆ makeBcString()

std::string CdrSolver::makeBcString ( const int  a_dir,
const Side::LoHiSide  a_side 
) const
protectedvirtual

Shortcut for making a boundary condition string.

Parameters
[in]a_dirDirection.
[in]a_sideCoordinate side.
Returns
Returns string of type m_className.bc.direction.side.

◆ nonConservativeDivergence()

void CdrSolver::nonConservativeDivergence ( EBAMRIVData a_nonConservativeDivergence,
const EBAMRCellData a_divG 
)
protectedvirtual

Compute the non-conservative divergence.

Parameters
[out]a_nonConservativeDivergenceNon-conservative divergence.
[in]a_divGConservative divergence.
Note
This will compute the non-conservative divergence as divNC(G) = sum(kappa*div(G))/sum(kappa) in the cut-cell neighborhood.

◆ operator=() [1/2]

CdrSolver& CdrSolver::operator= ( const CdrSolver &&  a_other)
delete

Disallowed move assignement operator.

Parameters
[in]a_otherOther solver

◆ operator=() [2/2]

CdrSolver& CdrSolver::operator= ( const CdrSolver a_other)
delete

Disallowed assignment operator.

Parameters
[in]a_otherOther solver

◆ preRegrid()

void CdrSolver::preRegrid ( const int  a_lbase,
const int  a_oldFinestLevel 
)
virtual

Perform pre-regrid operations.

Parameters
[in]a_lbaseCoarsest level that changed during regrid.
[in]a_oldFinestLevelFinest grid level before the regrid operation.
Note
This stores m_phi and m_source.

Reimplemented in CdrMultigrid.

◆ redistribute()

void CdrSolver::redistribute ( EBAMRCellData a_phi,
const EBAMRIVData a_delta 
) const
virtualnoexcept

Add data through redistribution into cell-centered holders.

Parameters
[in,out]a_phiTarget data
[in]a_deltaData to be redistributed

◆ registerOperators()

void CdrSolver::registerOperators ( )
virtual

Register operators for AMR operations.

Note
This includes operators for redistribution, flux registers, regridding, ghost cell interpolation, and conservative coarsening.

Reimplemented in CdrMultigrid.

◆ regrid()

void CdrSolver::regrid ( const int  a_lmin,
const int  a_oldFinestLevel,
const int  a_newFinestLevel 
)
virtual

Write checkpoint data into HDF5 file. @paramo[out] a_handle HDF5 file.

Parameters
[in]a_levelGrid level

Read checkpoint data from HDF5 file.

Parameters
[in]a_handleHDF5 handle.
[in]constint a_level Grid level

Regrid this solver.

Parameters
[in]a_lminCoarsest level where grids did not change.
[in]a_oldFinestLevelFinest AMR level before the regrid.
[in]a_newFinestLevelFinest AMR level after the regrid.

This linearly interpolates (with limiters) m_state to the new grids. Velocities are left undefined.

◆ resetDomainFlux()

void CdrSolver::resetDomainFlux ( EBAMRFluxData a_flux)
protectedvirtual

Set flux to zero on domain boundaries.

Parameters
[in]a_fluxFlux data holder – on domain edges this is modified so the flux is zero.

◆ setAmr()

void CdrSolver::setAmr ( const RefCountedPtr< AmrMesh > &  a_amr)
virtual

Set the amr object.

Parameters
[in]a_amrAmrMesh object.

◆ setComputationalGeometry()

void CdrSolver::setComputationalGeometry ( const RefCountedPtr< ComputationalGeometry > &  a_computationalGeometry)
virtual

Set computational geometry.

Parameters
[in]a_computationalGeometryThe computational geometry.

◆ setDefaultDomainBC()

void CdrSolver::setDefaultDomainBC ( )

This sets default boundary conditions (wall type).

Called in constructor for setting domain BCs.

Note
This sets the domain boundary condition to be a wall BC (no incoming/outgoing mass). This is the default behavior – if the user wants different BCs he can call setDomainBcType and setDomainBcFunction.

◆ setDiffusionCoefficient() [1/3]

void CdrSolver::setDiffusionCoefficient ( const EBAMRFluxData a_diffusionCoefficient,
const EBAMRIVData a_ebDiffusionCoefficient 
)
virtual

Data-based version of setting diffusion coefficients (which are stored on faces)

Parameters
[in]a_diffusionCoefficientFace-centered diffusion coefficient.
[in]a_ebDiffusionCoefficientEB-centered diffusion coefficient.
Note
Realms do not have to be the same.

◆ setDiffusionCoefficient() [2/3]

void CdrSolver::setDiffusionCoefficient ( const Real  a_diffusionCoefficient)
virtual

Set diffusion coefficient to be constant everywhere (they are stored on faces)

Parameters
[in]a_diffusionCoefficientGlobal diffusion coefficient.

This sets the diffusion coefficient to be constant on both faces and EB faces.

◆ setDiffusionCoefficient() [3/3]

void CdrSolver::setDiffusionCoefficient ( const std::function< Real(const RealVect a_position)> &  a_diffusionCoefficient)
virtual

Polymorphic way of setting diffusion coefficients every.

Parameters
[in]a_diffusionCoefficientDiffusion coefficient.
Note
This will evaluate the function in all grid cells and set the diffusion coefficient

◆ setDomainBcFunction()

void CdrSolver::setDomainBcFunction ( const CdrDomainBC::DomainSide  a_domainSide,
const CdrDomainBC::FluxFunction  a_fluxFunction 
)

Set domain bc function on particular domain side.

Parameters
[in]a_domainSideDomain side
[in]a_fluxFunctionFlux function
Note
This only associates a function with a domain side – to use it one must also call setDomainBcType(a_domainSide, BcType::Function)

◆ setDomainBcType()

void CdrSolver::setDomainBcType ( const CdrDomainBC::DomainSide  a_domainSide,
const CdrDomainBC::BcType  a_bcType 
)

Set domain bc type on domain side.

Parameters
[in]a_domainSideDomain side
[in]a_bcTypeBoundary condition type

◆ setEbFlux() [1/2]

void CdrSolver::setEbFlux ( const EBAMRIVData a_ebFlux)
virtual

Data-based version of setting the EB flux.

Parameters
[in]a_ebFluxFlux on EB centroids.
Note
a_ebFlux does not have to be the same realm as in the solver.

◆ setEbFlux() [2/2]

void CdrSolver::setEbFlux ( const Real  a_ebFlux)
virtual

Set the eb flux to a constant.

Parameters
[in]a_ebFluxFlux on EB centroids.

This sets a constant EB-flux on all EBs.

◆ setPhase()

void CdrSolver::setPhase ( phase::which_phase  a_phase)
virtual

Set phase.

Parameters
[in]a_phasePhase.

◆ setRealm()

void CdrSolver::setRealm ( const std::string  a_realm)
virtual

Set the realm for this solver.

Parameters
[in]a_realmRealm identifier

◆ setSource() [1/3]

void CdrSolver::setSource ( const EBAMRCellData a_source)
virtual

Data based version of setting source terms.

Parameters
[in]a_sourceSource term.

This sets m_source to be equal to a_source

◆ setSource() [2/3]

void CdrSolver::setSource ( const Real  a_source)
virtual

Set constant source terms everywhere.

Parameters
[in]a_sourceSource term.
Note
This sets m_source to a_source everywhere

◆ setSource() [3/3]

void CdrSolver::setSource ( const std::function< Real(const RealVect a_position)>  a_source)
virtual

Polymorphic way of setting source terms.

Parameters
[in]a_sourceSource term.
Note
This fills m_source with the evaluation of a_source

◆ setSpecies()

void CdrSolver::setSpecies ( const RefCountedPtr< CdrSpecies > &  a_species)
virtual

Set species.

Parameters
[in]a_speciesSpecies (where this solver gets its name, charge, initial conditions etc).

◆ setTime()

void CdrSolver::setTime ( const int  a_step,
const Real  a_time,
const Real  a_dt 
)
virtual

Set the time for this solver.

Parameters
[in]a_stepTime step number
[in]a_timeTime (in seconds)
[in]a_dtTime step increment
Note
This sets m_step=a_step, m_time=a_time, m_dt=a_dt

◆ setVelocity() [1/3]

void CdrSolver::setVelocity ( const EBAMRCellData a_velocity)
virtual

Set velocity from data holder.

Parameters
[in]a_velocityCell centered velocities.

This copies a_velo onto m_cellVelocity.

◆ setVelocity() [2/3]

void CdrSolver::setVelocity ( const RealVect  a_velocity)
virtual

Set constant velocity.

Parameters
[in]a_velocityVelocity

This sets the cell centered velocity to be a_velo in every grid cell.

◆ setVelocity() [3/3]

void CdrSolver::setVelocity ( const std::function< RealVect(const RealVect a_pos)> &  a_velocity)
virtual

Set velocity using a polymorphic function.

Parameters
[in]a_velocityVelocity to be set. You can use e.g. a lambda to set it.

◆ setVerbosity()

void CdrSolver::setVerbosity ( const int  a_verbosity)
virtual

Set verbosity.

Parameters
[in]a_verbosityVerbosity level.

◆ smoothHeavisideFaces()

void CdrSolver::smoothHeavisideFaces ( EBAMRFluxData a_facePhi,
const EBAMRCellData a_cellPhi 
)
protectedvirtual

Use Heaviside smoothing for computing face-centered states.

Parameters
[out]a_facePhiFace centered state
[in]a_cellPhiCell-centered states
Note
This is only intended to be used with FHD functionality.

◆ weightedUpwind()

void CdrSolver::weightedUpwind ( EBAMRCellData a_weightedUpwindPhi,
const int  a_pow 
)
virtual

Compute an upwind-weighted version of phi.

This is equal to sum(alpha * v^p * phi)/sum(alpha * v^p) where alpha is the face aperture and v is the face velocity. phi is taken as the upwind value.

Parameters
[out]a_weightedUpwindPhiWeighted upwind phi
[in]a_powWeighting factor

◆ writeData()

void CdrSolver::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
protectedvirtualnoexcept

Write data to output. Convenience function.

Parameters
[in,out]a_outputOutput data holder.
[in,out]a_icompStarting component where this solver begins writing the output.
[in]a_dataData to write.
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level
[in]a_interpToCentroidsIf true, a_data will be interpolated to cell centroids before writing to a_output.
[in]a_interpGhostIf true, interpolate ghost cells

◆ writePlotData()

void CdrSolver::writePlotData ( LevelData< EBCellFAB > &  a_output,
int &  a_icomp,
const std::string  a_outputRealm,
const int  a_level 
) const
virtualnoexcept

Write output data to a_output.

Parameters
[in,out]a_outputOutput data holder.
[in,out]a_icompStarting component where this solver begins writing the output.
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level
Note
This will write the plot data in this solver to a_output, starting on a_comp
This routine writes m_phi on centroids (and not cell-centers).

◆ writePlotFile()

void CdrSolver::writePlotFile ( )
virtual

Write plot file.

The name of the plot file is m_name.stepXXXXX.DIM.hdf5

Note
This calls writePlotData(...)

Member Data Documentation

◆ m_cellCenteredDiffusionCoefficient

EBAMRCellData CdrSolver::m_cellCenteredDiffusionCoefficient
protected

Diffusion coefficients on cell centers.

Note
Memory is only allocated if the solver is diffusive

◆ m_cellVelocity

EBAMRCellData CdrSolver::m_cellVelocity
protected

Cell-centered velocities.

Note
Memory is only allocated if the solver is mobile.

◆ m_className

std::string CdrSolver::m_className
protected

Class name.

This will be different for different implementations of CdrSolver

◆ m_domainFluxFunctions

std::map<CdrDomainBC::DomainSide, CdrDomainBC::FluxFunction> CdrSolver::m_domainFluxFunctions
protected

Domain flux functions.

Used when the user specifies that the domain flux is function-based

◆ m_ebCenteredDiffusionCoefficient

EBAMRIVData CdrSolver::m_ebCenteredDiffusionCoefficient
protected

Diffusion coefficients on EB faces.

Note
Memory is only allocated if the solver is diffusive

◆ m_faceCenteredDiffusionCoefficient

EBAMRFluxData CdrSolver::m_faceCenteredDiffusionCoefficient
protected

Diffusion coefficients on face centers.

Note
Memory is only allocated if the solver is diffusive

◆ m_faceStates

EBAMRFluxData CdrSolver::m_faceStates
protected

Holder for face centered states.

Note
Memory is only allocated if the solver is mobile.

◆ m_faceVelocity

EBAMRFluxData CdrSolver::m_faceVelocity
protected

Face-centered velocities (only normal components)

Note
Memory is only allocated if the solver is mobile.

◆ m_verbosity

int CdrSolver::m_verbosity
protected

Solver verbosity.

A high number gives lots of output


The documentation for this class was generated from the following files: