chombo-discharge
|
Base class for solving convection-diffusion-reaction equations. More...
#include <CD_CdrSolver.H>
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... | |
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. | |
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 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 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< 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. | |
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. | |
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).
CdrSolver::CdrSolver | ( | ) |
Default constructor.
|
delete |
Disallowed copy constructor.
[in] | a_other | Other solver |
|
delete |
Disallowed move constructor.
[in] | a_other | Other solver |
|
pure virtual |
Implicit diffusion Crank-Nicholson advance with source term.
[in,out] | a_newPhi | Solution at time t + dt |
[in] | a_oldPhi | Solution at time t |
[in] | a_source | Source term. |
[in] | a_dt | Time step |
Implemented in CdrMultigrid.
|
virtual |
Implicit diffusion Crank-Nicholson advance without source term.
[in,out] | a_newPhi | Solution at time t + dt |
[in] | a_oldPhi | Solution at time t |
[in] | a_dt | Time step |
|
pure virtual |
Implicit diffusion Euler advance with source term.
[in,out] | a_newPhi | Solution at time t + dt |
[in] | a_oldPhi | Solution at time t |
[in] | a_source | Source term. |
[in] | a_dt | Time step |
This solves the implicit diffusion equation equation a_newPhi - a_oldPhi = dt*Laplacian(a_newPhi) + dt*a_source.
Implemented in CdrMultigrid.
|
virtual |
Implicit diffusion Euler advance without source term.
[in,out] | a_newPhi | Solution at time t + dt |
[in] | a_oldPhi | Solution at time t |
[in] | a_dt | Time step |
|
protectedpure virtual |
Advection-only extrapolation to faces.
[out] | a_facePhi | Phi on faces |
[in] | a_phi | Phi on cell center |
[in] | a_extrapDt | Time centering/extrapolation (if the advective integrator can do it) |
Implemented in CdrMultigrid, CdrGodunov, and CdrCTU.
|
virtual |
Allocate internal storage.
Reimplemented in CdrMultigrid, and CdrGodunov.
|
virtual |
Average velocities to faces.
|
protectedvirtual |
Average cell-centered velocities to face centers.
[out] | a_faceVelocity | Face-centered velocities |
[in] | a_cellVelocity | Cell-centered velocities |
|
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.
|
protectedvirtual |
Compute the full advection-diffusion flux. This assumes that the solver is mobile and diffusive.
[out] | a_flux | Face-centered flux. |
[in] | a_cellStates | Cell-centered states. |
[in] | a_faceVelocities | Face-centered velocities. |
[in] | a_faceStates | Face-centered (advected) states. |
[in] | a_faceDiffCo | Face-centered diffusion coefficients. |
[in] | a_addDomainFlux | Add domain flux or not. |
|
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.
Reimplemented in CdrGodunov, and CdrCTU.
|
protectedvirtual |
Set up face-centered advection flux.
[out] | a_flux | Face-centered fluxes |
[in] | a_facePhi | Face-centered state |
[in] | a_faceVelocity | Face-centered velocities |
[in] | a_addDomainFlux | Add the domain flux |
|
protectedvirtual |
Set up face-centered advection flux on a grid level.
[out] | a_flux | Face-centered fluxes |
[in] | a_facePhi | Face-centered state |
[in] | a_faceVelocity | Face-centered velocities |
[in] | a_domainFlux | Domain flux |
[in] | a_level | Grid level |
|
virtual |
Compute the total charge in m_phi.
Calls computeMass and multiplies by charge
|
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.
|
protectedvirtual |
Compute the face-centered diffusion flux.
[out] | a_flux | Face-centered flux |
[in] | a_phi | Cell centered phi |
[in] | a_addDomainFlux | Add domain flux or not |
|
protectedvirtual |
Compute the face-centered diffusion flux.
[out] | a_flux | Face-centered flux |
[in] | a_phi | Cell centered phi |
[in] | a_lvl | Grid level |
|
pure virtual |
Compute div(D*grad(phi)) explicitly.
[out] | a_divF | Divergence term, i.e. finite volume approximation to Div(D*Grad(phi)). |
[in] | a_phi | Cell-centered state |
[in] | a_domainBc | Flag for setting domain boundary conditions |
[in] | a_conservativeOnly | If true, we compute div(D) = 1/dx*sum(fluxes), which does not involve redistribution. |
[in] | a_useEbFlux | If true, the embedded boundary flux will be injected and included in div(D) |
[in] | a_domainFlux | If true, the domain flux will injected and included in div(D) |
Implemented in CdrMultigrid.
|
protectedvirtual |
Compute conservative divergence on irregular cells (not kappa divided)
[out] | a_divG | Conservative divergence (not kappa divided) |
[out] | a_centroidFlux | Face centroid centered fluxes |
[in] | a_ebFlux | EB flux |
[in] | a_lvl | Grid level |
|
pure virtual |
Compute div(v*phi) explicitly.
[out] | a_divF | Divergence term, i.e. finite volume approximation to Div(v*phi), including redistribution magic. |
[in] | a_phi | Cell-centered state |
[in] | a_extrapDt | Extrapolation in time, i.e. shifting of div(F) towards e.g. half time step. Only affects the advective term. |
[in] | a_conservativeOnly | If true, we compute div(F)=1/dx*sum(fluxes), which does not involve redistribution. |
[in] | a_domainBc | How to set domain fluxes |
[in] | a_ebFlux | If true, the embedded boundary flux will injected be included in div(F) |
[in] | a_domainFlux | If true, the domain flux will injected and included in div(F) |
Implemented in CdrMultigrid.
|
virtual |
Compute div(G) where G is a general face-centered flux on face centers and EB centers. This can involve mass redistribution.
[in] | a_divG | div(G) or kappa*div(G). |
[in,out] | a_G | Vector field which contains face-centered fluxes on input. Contains face-centroid fluxes on output. |
[in] | a_ebFlux | Flux on the EB centroids |
[in] | a_conservativeOnly | If true, we compute div(G)=1/dx*sum(fluxes), which does not involve redistribution. |
|
pure virtual |
Compute div(J) explicitly, where J = nV - D*grad(n)
[out] | a_divJ | Divergence term, i.e. finite volume approximation to |
[in] | a_phi | Cell-centered state |
[in] | a_extrapDt | Extrapolation in time, i.e. shifting of div(J) towards e.g. half time step. Only affects the advective term. |
[in] | a_conservativeOnly | If true, we compute div(J) = 1/dx*sum(fluxes), which does not involve redistribution. |
[in] | a_ebFlux | If true, the embedded boundary flux will injected and included in div(J) |
[in] | a_domainFlux | If true, the domain flux will injected and included in div(J) |
Implemented in CdrMultigrid.
|
virtual |
Compute the "physical mass" in m_phi.
This conservatively coarsens the solution and computes on the coarsest level only
|
virtual |
Compute the "physical mass" in the input argument.
This conservatively coarsens the solution and computes on the coarsest level only
[in] | a_phi | Cell-centered data. |
[in] | a_kappaScale | Multiply by kappa in irregular cells. |
|
virtual |
Compute the largest possible source time step (for explicit methods.
[in] | a_max | Maximum value of m_phi |
[in] | a_tolerance | Tolerance |
This computes dt = phi/source, but only for cells where phi lies within a_tolerance*a_max (and source > 0)
|
protectedvirtual |
Compute conservative divergence from fluxes.
[out] | a_conservativeDivergence | Conservative divergence computed as div(F) using finite volumes. Not divided by kappa. |
[in] | a_flux | Face-centered fluxes. Includes domain fluxes. |
[in] | a_ebFlux | EB flux. |
|
protectedvirtual |
Compute the conservative divergence over regular cells.
[out] | a_divJ | Conservative divergence. Not kappa-divided. |
[in] | a_flux | Face-centered fluxes |
[in] | a_lvl | Grid level |
|
virtual |
Deallocate internal storage.
|
protectedvirtual |
Define stencils for doing face-centered to face-centroid-centered states.
|
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.
|
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.
[out] | a_ebFlux | Projected flux on the EB. |
|
protectedvirtual |
Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux.
[in,out] | a_flux | Flux to be modified. |
|
protectedvirtual |
Set domain in data holder. This sets the flux on the boundary to either zero or to m_domainFlux.
[in,out] | a_flux | Flux to be modified. |
[in] | a_level | Grid level |
|
protectedvirtual |
Gaussian noise field.
[out] | a_noise | Gaussian white nosie |
[out] | a_sigma | Standard deviation |
|
virtual |
Get the cell-centered velocity.
|
virtual |
Get the domain flux data holder.
|
virtual |
Get the EB-centered diffusion coefficient.
|
virtual |
Get the eb-centered velocities.
|
virtual |
Get the eb flux data holder.
|
virtual |
Get the face-centered diffusion coefficient.
|
virtual |
Get the face-centered velocities.
|
virtual |
Get solver name.
|
virtual |
Get number of output fields.
|
virtual |
Get the cell-centered phi.
|
virtual |
Get output plot names.
|
virtual |
Get the realm where this solver is registered.
|
virtual |
Get the source term.
|
virtual |
Compute a random gaussian white noise source term.
[out] | a_noiseSource | Source term |
[in] | a_cellPhi | Cell-centered states |
|
protectedvirtual |
Use the non-conservative divergence to make the conservative divergence hold the hybrid divergence.
[in,out] | a_hybridDivergence | On input, contains the conservative divergence (without kappa division). Contains hybrid divergence on output. |
[out] | a_massDifference | Mass difference between updating with hybrid divergence and true divergence |
[in] | a_nonConservativeDivergence | Non-conservative divergence. |
|
protectedvirtual |
Make the hybrid divergence. On the way in, a_hybridDivergence must hold the conservative divergence.
[in,out] | a_hybridDivergence | On input, contains the conservative divergence (without kappa division). Contains hybrid divergence on output. |
[out] | a_massDifference | Mass difference between updating with hybrid divergence and true divergence |
[in] | a_nonConservativeDivergence | Non-conservative divergence. |
[in] | a_lvl | Grid level. |
|
virtual |
Fill m_phi state with initial data from m_species.
|
protectedvirtual |
Fill initial data from a distribution function.
|
protectedvirtual |
Fill initial data from particles.
|
protectedvirtual |
Interpolate flux to centroids.
[in,out] | a_flux | On input, contains centered fluxes. Output contains centroid fluxes |
This interpolates face fluxes to centroids.
|
protectedvirtual |
Interpolate flux to centroids.
[in,out] | a_flux | On input, contains centered fluxes. Output contains centroid fluxes |
[in] | a_lvl | Grid level |
This interpolates face fluxes to centroids.
|
protectedvirtual |
Shortcut for making a boundary condition string.
[in] | a_dir | Direction. |
[in] | a_side | Coordinate side. |
|
protectedvirtual |
Compute the non-conservative divergence.
[out] | a_nonConservativeDivergence | Non-conservative divergence. |
[in] | a_divG | Conservative divergence. |
Disallowed move assignement operator.
[in] | a_other | Other solver |
Disallowed assignment operator.
[in] | a_other | Other solver |
|
virtual |
Perform pre-regrid operations.
[in] | a_lbase | Coarsest level that changed during regrid. |
[in] | a_oldFinestLevel | Finest grid level before the regrid operation. |
Reimplemented in CdrMultigrid.
|
virtualnoexcept |
Add data through redistribution into cell-centered holders.
[in,out] | a_phi | Target data |
[in] | a_delta | Data to be redistributed |
|
virtual |
Register operators for AMR operations.
Reimplemented in CdrMultigrid.
|
virtual |
Write checkpoint data into HDF5 file. @paramo[out] a_handle HDF5 file.
[in] | a_level | Grid level |
Read checkpoint data from HDF5 file.
[in] | a_handle | HDF5 handle. |
[in] | const | int a_level Grid level |
Regrid this solver.
[in] | a_lmin | Coarsest level where grids did not change. |
[in] | a_oldFinestLevel | Finest AMR level before the regrid. |
[in] | a_newFinestLevel | Finest AMR level after the regrid. |
This linearly interpolates (with limiters) m_state to the new grids. Velocities are left undefined.
|
protectedvirtual |
Set flux to zero on domain boundaries.
[in] | a_flux | Flux data holder – on domain edges this is modified so the flux is zero. |
|
virtual |
Set the amr object.
[in] | a_amr | AmrMesh object. |
|
virtual |
Set computational geometry.
[in] | a_computationalGeometry | The computational geometry. |
void CdrSolver::setDefaultDomainBC | ( | ) |
This sets default boundary conditions (wall type).
Called in constructor for setting domain BCs.
|
virtual |
Data-based version of setting diffusion coefficients (which are stored on faces)
[in] | a_diffusionCoefficient | Face-centered diffusion coefficient. |
[in] | a_ebDiffusionCoefficient | EB-centered diffusion coefficient. |
|
virtual |
Set diffusion coefficient to be constant everywhere (they are stored on faces)
[in] | a_diffusionCoefficient | Global diffusion coefficient. |
This sets the diffusion coefficient to be constant on both faces and EB faces.
|
virtual |
Polymorphic way of setting diffusion coefficients every.
[in] | a_diffusionCoefficient | Diffusion coefficient. |
void CdrSolver::setDomainBcFunction | ( | const CdrDomainBC::DomainSide | a_domainSide, |
const CdrDomainBC::FluxFunction | a_fluxFunction | ||
) |
Set domain bc function on particular domain side.
[in] | a_domainSide | Domain side |
[in] | a_fluxFunction | Flux function |
void CdrSolver::setDomainBcType | ( | const CdrDomainBC::DomainSide | a_domainSide, |
const CdrDomainBC::BcType | a_bcType | ||
) |
Set domain bc type on domain side.
[in] | a_domainSide | Domain side |
[in] | a_bcType | Boundary condition type |
|
virtual |
Data-based version of setting the EB flux.
[in] | a_ebFlux | Flux on EB centroids. |
|
virtual |
Set the eb flux to a constant.
[in] | a_ebFlux | Flux on EB centroids. |
This sets a constant EB-flux on all EBs.
|
virtual |
Set phase.
[in] | a_phase | Phase. |
|
virtual |
Set the realm for this solver.
[in] | a_realm | Realm identifier |
|
virtual |
Data based version of setting source terms.
[in] | a_source | Source term. |
This sets m_source to be equal to a_source
|
virtual |
Set constant source terms everywhere.
[in] | a_source | Source term. |
|
virtual |
Polymorphic way of setting source terms.
[in] | a_source | Source term. |
|
virtual |
Set species.
[in] | a_species | Species (where this solver gets its name, charge, initial conditions etc). |
|
virtual |
Set the time for this solver.
[in] | a_step | Time step number |
[in] | a_time | Time (in seconds) |
[in] | a_dt | Time step increment |
|
virtual |
Set velocity from data holder.
[in] | a_velocity | Cell centered velocities. |
This copies a_velo onto m_cellVelocity.
|
virtual |
Set constant velocity.
[in] | a_velocity | Velocity |
This sets the cell centered velocity to be a_velo in every grid cell.
|
virtual |
Set velocity using a polymorphic function.
[in] | a_velocity | Velocity to be set. You can use e.g. a lambda to set it. |
|
virtual |
Set verbosity.
[in] | a_verbosity | Verbosity level. |
|
protectedvirtual |
Use Heaviside smoothing for computing face-centered states.
[out] | a_facePhi | Face centered state |
[in] | a_cellPhi | Cell-centered states |
|
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.
[out] | a_weightedUpwindPhi | Weighted upwind phi |
[in] | a_pow | Weighting factor |
|
protectedvirtualnoexcept |
Write data to output. Convenience function.
[in,out] | a_output | Output data holder. |
[in,out] | a_icomp | Starting component where this solver begins writing the output. |
[in] | a_data | Data to write. |
[in] | a_outputRealm | Realm where a_output belongs |
[in] | a_level | Grid level |
[in] | a_interpToCentroids | If true, a_data will be interpolated to cell centroids before writing to a_output. |
[in] | a_interpGhost | If true, interpolate ghost cells |
|
virtualnoexcept |
Write output data to a_output.
[in,out] | a_output | Output data holder. |
[in,out] | a_icomp | Starting component where this solver begins writing the output. |
[in] | a_outputRealm | Realm where a_output belongs |
[in] | a_level | Grid level |
|
virtual |
Write plot file.
The name of the plot file is m_name.stepXXXXX.DIM.hdf5
|
protected |
Diffusion coefficients on cell centers.
|
protected |
Cell-centered velocities.
|
protected |
Class name.
This will be different for different implementations of CdrSolver
|
protected |
Domain flux functions.
Used when the user specifies that the domain flux is function-based
|
protected |
Diffusion coefficients on EB faces.
|
protected |
Diffusion coefficients on face centers.
|
protected |
Holder for face centered states.
|
protected |
Face-centered velocities (only normal components)
|
protected |
Solver verbosity.
A high number gives lots of output