Class Maestro

class Maestro : public amrex::AmrCore

Public Types

enum plotfile_flag

flag for writing plotfiles

Values:

enumerator plotInitData
enumerator plotInitProj
enumerator plotDivuIter

Public Functions

Maestro()

constructor

~Maestro() override

destructor

Maestro(Maestro &other) = delete
Maestro(Maestro &&other) = delete
Maestro &operator=(const Maestro &other) = delete
Maestro &operator=(Maestro &&other) = delete
void Setup()

Setup the simulation.

  • read in C++/F90 parameters

  • define global C++/F90 variables and initialize network

  • set up boundary conditions

  • initialize base state geometry parameters

  • set istep, t_new, t_old

  • allocate MultiFabs and base state arrays

void Init()

Initialize the simulation.

  • initialize multifab and base state data

  • perform initial projection

  • perform divu iters

  • perform initial (pressure) iterations

void Evolve()

advance solution to final time

void AdvanceTimeStep(bool is_initIter)

Advance solution at all levels for a single time step. This uses the old temporal integration scheme.

Parameters:

is_initIter – is it the initial iteration?

void AdvanceTimeStepIrreg(bool is_initIter)

Advance solution for a single time step with irregular base state spacing

Parameters:

is_initIter – is it the initial iteration?

void AdvanceTimeStepAverage(bool is_initIter)

Advance solution for a single time step with regular base state spacing and new time-stepping scheme

Parameters:

is_initIter – is it the initial iteration?

void AdvectBaseDens(BaseState<amrex::Real> &rho0_predicted_edge)
void AdvectBaseDensPlanar(BaseState<amrex::Real> &rho0_predicted_edge_state)
void AdvectBaseDensSphr(BaseState<amrex::Real> &rho0_predicted_edge_state)
void AdvectBaseEnthalpy(BaseState<amrex::Real> &rhoh0_predicted_edge)
void AdvectBaseEnthalpyPlanar(BaseState<amrex::Real> &rhoh0_predicted_edge_state)
void AdvectBaseEnthalpySphr(BaseState<amrex::Real> &rhoh0_predicted_edge_state)
void UpdateSpecies(const BaseState<Real> &rho0, const BaseState<Real> &rho0_predicted_edge, const BaseState<Real> &rhoX0_old, BaseState<Real> &rhoX0_new)
void AdvancePremac(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const amrex::Vector<amrex::MultiFab> &w0_force_cart)

Compute unprojected mac velocities

Parameters:
  • umac – MAC velocity

  • w0mac – MAC base-state velocity

  • w0_force – base-state-velocity force

  • w0_force_cart – base-state-velocity force on cartesian grid

void MakeRhoXFlux(const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sflux, amrex::Vector<amrex::MultiFab> &etarhoflux, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_edge_old_state, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &r0mac_old, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &rho0_edge_new_state, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &r0mac_new, const BaseState<amrex::Real> &rho0_predicted_edge_state, int start_comp, int num_comp)

Calculate rhoX flux

Takes the predicted edges states of the scalars and the MAC velocities and computes the flux through the interfaces.

The construction of the fluxes depends on what form the incoming edge states take. This depends on species_pred_type:

  • predict_rhoprime_and_X: We have rho’ and X, and need a edge-centered base state to make the final fluxes

  • predict_rhoX: We use the (rho X) edge state directly to compute the fluxes. No base state input needed.

  • predict_rho_and_X: The fluxes are computed from the product of the rho and X edge states, again, no base state input needed.

Parameters:
  • state – cell-centered scalars

  • sflux – scalar flux

  • etarhofluxeta_rho flux

  • sedge – edge state of scalars

  • umac – MAC velocity of full state

  • r0_old – old base-state density

  • r0_edge_old – old base-state density on cell-edges

  • r0mac_old – old MAC-projected base-state density

  • r0_new – new base-state density

  • r0_edge_new – new base-state density on cell-edges

  • r0mac_new – new MAC-projected base-state density

  • r0_predicted_edge – new base-state density on cell edges

  • start_comp – index of component of state to begin with

  • num_comp – number of components to perform calculation for

void MakeRhoHFlux(const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sflux, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_edge_old, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &r0mac_old, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &rho0_edge_new, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &r0mac_new, const BaseState<amrex::Real> &rhoh0_old_in, const BaseState<amrex::Real> &rhoh0_edge_old, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &rh0mac_old, const BaseState<amrex::Real> &rhoh0_new_in, const BaseState<amrex::Real> &rhoh0_edge_new, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &rh0mac_new, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &h0mac_old, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &h0mac_new)

Calculate rhoh flux

Takes the predicted edges states of the scalars and the MAC velocities and computes the flux through the interfaces.

Parameters:
  • state – cell-centered scalars

  • sflux – scalar flux

  • sedge – edge state of scalars

  • umac, w0mac – MAC velocity of full and base-state velocity

  • r0_old – old base-state density

  • r0_edge_old – old base-state density on cell-edges

  • r0mac_old – old MAC-projected base-state density

  • r0_new – new base-state density

  • r0_edge_new – new base-state density on cell-edges

  • r0mac_new – new MAC-projected base-state density

  • rh0_old – old base-state conserved enthalpy

  • rh0_edge_old – old base-state conserved enthalpy on cell-edges

  • rh0mac_old – old MAC-projected base-state conserved enthalpy

  • rh0_new – new base-state conserved enthalpy

  • rh0_edge_new – new base-state conserved enthalpy on cell-edges

  • rh0mac_new – new MAC-projected base-state conserved enthalpy

  • h0mac_old, h0mac_new – base-state primitive enthalpy

void UpdateScal(const amrex::Vector<amrex::MultiFab> &stateold, amrex::Vector<amrex::MultiFab> &statenew, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sflux, const amrex::Vector<amrex::MultiFab> &force, int start_comp, int num_comp, const amrex::Vector<amrex::MultiFab> &p0_cart)

Given scalar fluxes, update scalars

Parameters:
  • stateold – cell-centered scalars

  • statenew – new scalar flux

  • sflux – scalar fluxes

  • force – velocity force

  • start_scomp – index of component of state to begin with

  • num_comp – number of components to perform calculation for

  • p0_cart – base state pressure on cartesian grid

void UpdateVel(const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &uedge, const amrex::Vector<amrex::MultiFab> &force, const amrex::Vector<amrex::MultiFab> &sponge, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac)

Update velocity

Parameters:
  • umac – MAC velocity

  • uedge – edge-based velocity

  • force – velocity force

  • sponge

  • w0mac – base state MAC velocity

void Average(const amrex::Vector<amrex::MultiFab> &phi, BaseState<amrex::Real> &phibar, int comp)

Compute the radial average of a quantity

Parameters:
  • mf – MultiFab containing quantity to be averaged

  • phibar – Averaged quantity

  • comp – Index of component of mf to average

void InitBaseState(BaseState<amrex::Real> &rho0, BaseState<amrex::Real> &rhoh0, BaseState<amrex::Real> &p0, const int lev)
void InitBaseStateMapSphr(const int lev, const amrex::MFIter &mfi, const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM> dx_fine, const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM> dx_lev)
void ComputeCutoffCoords(const BaseState<amrex::Real> &rho0_state) const
void RestrictBase(BaseState<Real> &s0, const bool is_cell_centered) const
void RestrictBase(const BaseStateArray<Real> &s0, const bool is_cell_centered) const
void FillGhostBase(BaseState<Real> &s0, const bool is_cell_centered) const
void FillGhostBase(const BaseStateArray<Real> &s0, const bool is_cell_centered) const
void CelltoEdge(const BaseState<amrex::Real> &s0_cell_s, BaseState<amrex::Real> &s0_edge_s) const
void WriteCheckPoint(int step)

Write a checkpoint at timestep step

int ReadCheckPoint()
void GotoNextLine(std::istream &is)
void PutInPertForm(amrex::Vector<amrex::MultiFab> &scal, const BaseState<Real> &s0, int comp, int bccomp, const amrex::Vector<amrex::BCRec> &bcs, bool flag)

If flag, subtract the base state, returning the perturbed quantity. Otherwise, add the base state, returning the full quantity. This version iterates over all levels.

Parameters:
  • scal – full/perturbed scalar quantity to subtract/add base state to

  • s0 – base state scalar

  • comp – component of scal to perform calculation on

  • bccomp – component of bcs_in to use to enforce boundary conditions

  • bcs_in – boundary conditions

  • flag – determines whether base state is subtracted (true) or added (false)

void PutInPertForm(int level, amrex::Vector<amrex::MultiFab> &scal, const BaseState<Real> &s0, int comp, int bccomp, const amrex::Vector<amrex::BCRec> &bcs, bool flag)

If flag, subtract the base state, returning the perturbed quantity. Otherwise, add the base state, returning the full quantity. This version operates only on a single level.

Parameters:
  • level – level to perform calculation on

  • scal – full/perturbed scalar quantity to subtract/add base state to

  • s0 – base state scalar

  • comp – component of scal to perform calculation on

  • bccomp – component of bcs_in to use to enforce boundary conditions

  • bcs_in – boundary conditions

  • flag – determines whether base state is subtracted (true) or added (false)

void ConvertRhoXToX(amrex::Vector<amrex::MultiFab> &scal, bool flag)

If flag, returns species mass fraction X given the conserved variable rhoX. Otherwise, performs inverse operation

void ConvertRhoHToH(amrex::Vector<amrex::MultiFab> &scal, bool flag)

If flag, return enthalpy h given the conserved variable rhoh. Otherwise, performs inverse operation

void PrintBase(const RealVector &base, const bool is_cell_centered = true)

Print out the contents of the base state.

void PrintBase(const BaseState<amrex::Real> &base, const bool is_cell_centered = true)
void PrintMF(const amrex::Vector<amrex::MultiFab> &MF)

Print out the contents of a Vector of MultiFabs.

void PrintBA(const amrex::Vector<amrex::BoxArray> &ba)

Print out grid data from a Vector of BoxArrays.

void PrintEdge(const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &EDGE, int dir)

Print out the contents of a Vector of edge-based MultiFabs.

void WriteMF(const amrex::Vector<amrex::MultiFab> &mf, const std::string &name)

Utility to write out a multilevel multifab to a plotfile.

void DensityAdvance(int which_step, amrex::Vector<amrex::MultiFab> &scalold, amrex::Vector<amrex::MultiFab> &scalnew, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sflux, amrex::Vector<amrex::MultiFab> &scal_force, amrex::Vector<amrex::MultiFab> &etarhoflux, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const BaseState<amrex::Real> &rho0_predicted_edge)

Advance the density

Parameters:
  • which_step – Is this the predictor (1) or corrector (2) step?

  • scalold – old cell-centered scalars

  • scalnew – new scalars

  • sedge – edge-based scalars

  • sflux – scalar fluxes

  • scal_force – scalar force

  • etarhofluxeta_rho flux

  • umac – MAC velocity

  • w0mac – MAC base state velocity

  • rho0_predicted_edge – base state density predicted to cell edges

void WriteDiagFile(int &index)

Put together an array of multifabs for writing.

void DiagFile(const int step, const amrex::Real t_in, const BaseState<amrex::Real> &rho0_in, const BaseState<amrex::Real> &p0_in, const amrex::Vector<amrex::MultiFab> &u_in, const amrex::Vector<amrex::MultiFab> &s_in, int &index)

Write plotfile to disk.

void EstDt()

Compute the time step.

void FirstDt()

Compute initial time step.

void EstDt_Divu(BaseState<amrex::Real> &gp0, const BaseState<amrex::Real> &p0, const BaseState<amrex::Real> &gamma1bar) const
void EnforceHSE(const BaseState<amrex::Real> &rho0_s, BaseState<amrex::Real> &p0_s, const BaseState<amrex::Real> &grav_cell_s)
void EnthalpyAdvance(int which_step, amrex::Vector<amrex::MultiFab> &scalold, amrex::Vector<amrex::MultiFab> &scalnew, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sflux, amrex::Vector<amrex::MultiFab> &scal_force, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const amrex::Vector<amrex::MultiFab> &thermal)

Advance the enthalpy

Parameters:
  • which_step – Is this the predictor (1) or corrector (2) step?

  • scalold – old cell-centered scalars

  • scalnew – new scalars

  • sedge – edge-based scalars

  • sflux – scalar fluxes

  • scal_force – scalar force

  • umac – MAC velocity

  • w0mac – MAC base state velocity

  • thermal – thermal term

void FillPatch(amrex::Real time, amrex::Vector<amrex::MultiFab> &mf, amrex::Vector<amrex::MultiFab> &mf_old, amrex::Vector<amrex::MultiFab> &mf_new, int srccomp, int destcomp, int ncomp, int startbccomp, const amrex::Vector<amrex::BCRec> &bcs_in, int variable_type = 0)

Call FillPatch for all levels.

void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf, amrex::Vector<amrex::MultiFab> &mf_old, amrex::Vector<amrex::MultiFab> &mf_new, int srccomp, int destcomp, int ncomp, int startbccomp, const amrex::Vector<amrex::BCRec> &bcs_in, int variable_type = 0)

Compute a new multifab by coping in phi from valid region and filling ghost cells

  • works for single level and 2-level cases (fill fine grid ghost by interpolating from coarse)

  • srccomp is the source component

  • destcomp is the destination component AND the bc component

void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab &mf, amrex::Vector<amrex::MultiFab> &mf_old, amrex::Vector<amrex::MultiFab> &mf_new, int srccomp, int destcomp, int ncomp, const amrex::Vector<amrex::BCRec> &bcs, int variable_type = 0)

Fill an entire multifab by interpolating from the coarser level

  • this comes into play when a new level of refinement appears

  • srccomp is the source component

  • destcomp is the destination component AND the bc component

void GetData(int lev, amrex::Real time, amrex::Vector<amrex::MultiFab*> &mf, amrex::Vector<amrex::Real> &mftime, amrex::Vector<amrex::MultiFab> &mf_old, amrex::Vector<amrex::MultiFab> &mf_new)

Utility to copy in data from mf_old and/or mf_new into mf

  • if time=t_old we copy mf_old intomf

  • if time=t_new we copy mf_new into mf

  • otherwise copy in both mf_old and mf_new into mf and the fillpatch routines know to interpolate in time.

void AverageDown(amrex::Vector<amrex::MultiFab> &mf, int comp, int ncomp)

Set covered coarse cells to be the average of overlying fine cells

Parameters:
  • mf – MultiFab to average

  • comp – Index of first component to average

  • ncomp – Number of components to average

void AverageDownFaces(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &edge)

Set covered faces to be the average of overlying fine faces.

void FillUmacGhost(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac_in, int level = -1)

Fill in ONE ghost cell for all components of a face-centered (MAC) velocity field behind physical boundaries. Does not modify the velocities on the boundary

void FillPatchUedge(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &uedge)

Fill in all ghost cells for an edge-based MAC velocity field.

void Put1dArrayOnCart(const BaseState<amrex::Real> &s0, amrex::Vector<amrex::MultiFab> &s0_cart, const bool is_input_edge_centered, const bool is_output_a_vector, const amrex::Vector<amrex::BCRec> &bcs = amrex::Vector<amrex::BCRec>(), const int sbccomp = 0, const int variable_type = 0)

Maps 1d arrays onto multi-D cartesian MultiFabs

Parameters:
  • s0 – 1d base state

  • s0_cart – base state mapped to multi-d cartesian MultiFab

  • is_input_edge_centered – is the input edge-centered?

  • is_output_a_vector – is the output a vector?

  • bcs – boundary conditions

  • sbccomp – start boundary conditions component

void Put1dArrayOnCart(const int level, const BaseState<amrex::Real> &s0, amrex::MultiFab &s0_cart, const bool is_input_edge_centered, const bool is_output_a_vector, const amrex::Vector<amrex::BCRec> &bcs = amrex::Vector<amrex::BCRec>(), const int sbccomp = 0)

Maps 1d arrays onto multi-D cartesian MultiFabs

Parameters:
  • level – AMR level to perform calculation on

  • s0 – 1d base state

  • s0_cart – base state mapped to multi-d cartesian MultiFab

  • is_input_edge_centered – is the input edge-centered?

  • is_output_a_vector – is the output a vector?

  • bcs – boundary conditions

  • sbccomp – start boundary conditions component

void Put1dArrayOnCart(const int level, const BaseState<amrex::Real> &s0, amrex::Vector<amrex::MultiFab> &s0_cart, const bool is_input_edge_centered, const bool is_output_a_vector, const amrex::Vector<amrex::BCRec> &bcs = amrex::Vector<amrex::BCRec>(), const int sbccomp = 0)
void Addw0(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &u_edge, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const amrex::Real &mult)

Add (mult times) the MAC-projected base state velocity to the edge-based velocity uedge

Parameters:
  • uedge – edge based velocity

  • w0mac – MAC base state velocity

  • mult – multiplication factor

void MakeW0mac(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac)

MAC-project the base state velocity.

void MakeS0mac(const BaseState<amrex::Real> &s0, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &s0mac)

MAC-project the base state scalar s0

void MakeNormal()

Create the unit normal across the grids.

void PutDataOnFaces(const amrex::Vector<amrex::MultiFab> &s_cc, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &face, const bool harmonic_avg)

Put the cell-centered data s_cc on faces by averaging adjacent cells.

void MakeCCtoRadii()
void MakeVelForce(amrex::Vector<amrex::MultiFab> &vel_force_cart, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &uedge_in, const amrex::Vector<amrex::MultiFab> &rho, const BaseState<amrex::Real> &rho0, const BaseState<amrex::Real> &grav_cell, const amrex::Vector<amrex::MultiFab> &w0_force_cart, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const bool is_final_update, int do_add_utilde_force)

Calculate the velocity force term.

void ModifyScalForce(amrex::Vector<amrex::MultiFab> &scal_force, const amrex::Vector<amrex::MultiFab> &state, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac_in, const BaseState<amrex::Real> &s0_edge, const amrex::Vector<amrex::MultiFab> &s0_cart, int comp, const amrex::Vector<amrex::BCRec> &bcs, int fullform)
void MakeRhoHForce(amrex::Vector<amrex::MultiFab> &scal_force, const int is_prediction, const amrex::Vector<amrex::MultiFab> &thermal, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac_cart, const int add_thermal, const int &which_step)

Calculate the conserved enthalpy force term.

void MakeTempForce(amrex::Vector<amrex::MultiFab> &temp_force, const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &thermal, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac_in)
void MakeGamma1bar(const amrex::Vector<amrex::MultiFab> &scal, BaseState<amrex::Real> &gamma1bar, const BaseState<amrex::Real> &p0)

Calculate the horizontal average of \(\Gamma_1\).

void InitData()

fill in multifab and base state data

void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override

During initialization of a simulation, Maestro::InitData() calls AmrCore::InitFromScratch(), which calls a MakeNewGrids() function that repeatedly calls this function to create finer levels. This function creates a new fine level that did not exist before by interpolating from the coarser level overrides the pure virtual function in AmrCore

void InitProj()

Performs the initial projection.

void DivuIter(int istep_divu_iter)

Performs the divu iteration.

void InitIter()

Performs the initial iteration to initialize gradpi

void InitLevelData(const int lev, const amrex::Real time, const amrex::MFIter &mfi, const amrex::Array4<amrex::Real> scal, const amrex::Array4<amrex::Real> vel)
void InitLevelDataSphr(const int lev, const amrex::Real time, amrex::MultiFab &scal, amrex::MultiFab &vel)
void SetInletBCs()
void MacProj(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, amrex::Vector<amrex::MultiFab> &macphi, const amrex::Vector<amrex::MultiFab> &macrhs, const BaseState<amrex::Real> &beta0, const bool is_predictor)

Do the MAC projection

Parameters:
  • umac – enters with face-centered, time-centered Utilde^* and should leave with Utilde

  • macphi – is the solution to the elliptic solve and enters as either zero, or the solution to the predictor MAC projection

  • macrhs – enters as beta0*(S-Sbar)

  • beta0 – is a 1d cell-centered array

void MultFacesByBeta0(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &edge, const BaseState<amrex::Real> &beta0_s, const BaseState<amrex::Real> &beta0_edge_s, const int &mult_or_div)

Multiply (or divide) face-centered data by beta0

Parameters:
  • edge – face-centered data

  • beta0 – cell-centered \(\beta_0\)

  • beta0_edge – face-centered \(\beta_0\)

  • mult_or_div – do we multiply or divide?

void ComputeMACSolverRHS(amrex::Vector<amrex::MultiFab> &solverrhs, const amrex::Vector<amrex::MultiFab> &macrhs, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac)

Compute the RHS for the solve, RHS = macrhs - div(beta0*umac)

Parameters:
  • solverrhs – RHS for the solve

  • macrhsmacrhs term

  • umac – MAC velocity

void AvgFaceBcoeffsInv(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &facebcoef, const amrex::Vector<amrex::MultiFab> &rhocc)

Average bcoefs at faces using inverse of rho

Parameters:
  • facebcoef – face-centered bcoefs

  • rhocc – cell-centered density

void SetMacSolverBCs(amrex::MLABecLaplacian &mlabec)

Set boundaries for LABecLaplacian to solve -div(B grad) phi = RHS

void MakeBeta0(BaseState<amrex::Real> &beta0_s, const BaseState<amrex::Real> &rho0_s, const BaseState<amrex::Real> &p0_s, const BaseState<amrex::Real> &gamma1bar_s, const BaseState<amrex::Real> &grav_cell_s) const
void MakeEdgeScal(amrex::Vector<amrex::MultiFab> &state, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, amrex::Vector<amrex::MultiFab> &force, const bool is_vel, const amrex::Vector<amrex::BCRec> &bcs, int nbccomp, int start_scomp, int start_bccomp, int num_comp, const bool is_conservative)

Calculate scalars on cell edges given cell-centered state constructs the edge state of a scalar, using a second-order Taylor expansion in space (through dx/2) and time (though dt/2) (if ppm_type = 0) or using PPM (for ppm_type = 1,2).

We use only MAC-projected edge velocities in this prediction.

We are computing all edge states for each variable. This is what is done for the final updates of the state variables and velocity. For velocity, we should set is_vel = true

Parameters:
  • state – cell-centered scalars

  • sedge – edge state of scalars

  • umac – MAC velocity

  • force – velocity force

  • is_vel – set to true if state is a velocity

  • bcs – boundary conditions

  • nbccomp – number of components of bcs

  • start_scomp – index of component of state to begin with

  • start_bccomp – index of component of bcs to begin with

  • num_comp – number of components to perform calculation for

  • is_conservative – are these conserved quantities?

void MakeEdgeScalPredictor(const amrex::MFIter &mfi, amrex::Array4<amrex::Real> const slx, amrex::Array4<amrex::Real> const srx, amrex::Array4<amrex::Real> const sly, amrex::Array4<amrex::Real> const sry, amrex::Array4<amrex::Real> const scal, amrex::Array4<amrex::Real> const Ip, amrex::Array4<amrex::Real> const Im, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const simhx, amrex::Array4<amrex::Real> const simhy, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, int comp, int bccomp, bool is_vel) const
void MakeEdgeScalEdges(const amrex::MFIter &mfi, amrex::Array4<amrex::Real> const slx, amrex::Array4<amrex::Real> const srx, amrex::Array4<amrex::Real> const sly, amrex::Array4<amrex::Real> const sry, amrex::Array4<amrex::Real> const scal, amrex::Array4<amrex::Real> const sedgex, amrex::Array4<amrex::Real> const sedgey, amrex::Array4<amrex::Real> const force, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const Ipf, amrex::Array4<amrex::Real> const Imf, amrex::Array4<amrex::Real> const simhx, amrex::Array4<amrex::Real> const simhy, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, int comp, int bccomp, const bool is_vel, const bool is_conservative) const
void MakeDivU(const amrex::Box &bx, amrex::Array4<amrex::Real> const divu, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const wmac, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx) const
void MakeEdgeScalPredictor(const amrex::MFIter &mfi, amrex::Array4<amrex::Real> const slx, amrex::Array4<amrex::Real> const srx, amrex::Array4<amrex::Real> const sly, amrex::Array4<amrex::Real> const sry, amrex::Array4<amrex::Real> const slz, amrex::Array4<amrex::Real> const srz, amrex::Array4<amrex::Real> const scal, amrex::Array4<amrex::Real> const Ip, amrex::Array4<amrex::Real> const Im, amrex::Array4<amrex::Real> const slopez, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const wmac, amrex::Array4<amrex::Real> const simhx, amrex::Array4<amrex::Real> const simhy, amrex::Array4<amrex::Real> const simhz, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, int comp, int bccomp, const bool is_vel) const
void MakeEdgeScalTransverse(const amrex::MFIter &mfi, amrex::Array4<amrex::Real> const slx, amrex::Array4<amrex::Real> const srx, amrex::Array4<amrex::Real> const sly, amrex::Array4<amrex::Real> const sry, amrex::Array4<amrex::Real> const slz, amrex::Array4<amrex::Real> const srz, amrex::Array4<amrex::Real> const scal, amrex::Array4<amrex::Real> const divu, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const wmac, amrex::Array4<amrex::Real> const simhx, amrex::Array4<amrex::Real> const simhy, amrex::Array4<amrex::Real> const simhz, amrex::Array4<amrex::Real> const simhxy, amrex::Array4<amrex::Real> const simhxz, amrex::Array4<amrex::Real> const simhyx, amrex::Array4<amrex::Real> const simhyz, amrex::Array4<amrex::Real> const simhzx, amrex::Array4<amrex::Real> const simhzy, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, int comp, int bccomp, const bool is_vel, const bool is_conservative) const
void MakeEdgeScalEdges(const amrex::MFIter &mfi, amrex::Array4<amrex::Real> const slx, amrex::Array4<amrex::Real> const srx, amrex::Array4<amrex::Real> const sly, amrex::Array4<amrex::Real> const sry, amrex::Array4<amrex::Real> const slz, amrex::Array4<amrex::Real> const srz, amrex::Array4<amrex::Real> const scal, amrex::Array4<amrex::Real> const sedgex, amrex::Array4<amrex::Real> const sedgey, amrex::Array4<amrex::Real> const sedgez, amrex::Array4<amrex::Real> const force, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const wmac, amrex::Array4<amrex::Real> const Ipf, amrex::Array4<amrex::Real> const Imf, amrex::Array4<amrex::Real> const simhxy, amrex::Array4<amrex::Real> const simhxz, amrex::Array4<amrex::Real> const simhyx, amrex::Array4<amrex::Real> const simhyz, amrex::Array4<amrex::Real> const simhzx, amrex::Array4<amrex::Real> const simhzy, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, int comp, int bccomp, const bool is_vel, const bool is_conservative) const
void MakeEdgeState1d(BaseState<amrex::Real> &s, BaseState<amrex::Real> &sedge, BaseState<amrex::Real> &force)
void MakeEdgeState1dSphr(BaseState<amrex::Real> &s_state, BaseState<amrex::Real> &sedge_state, BaseState<amrex::Real> &force_state) const
void MakeEdgeState1dPlanar(BaseState<amrex::Real> &s_state, BaseState<amrex::Real> &sedge_state, BaseState<amrex::Real> &force_state) const
void MakeEtarho(const amrex::Vector<amrex::MultiFab> &etarho_flux)

Compute eta_rho at edge- and cell-centers

Parameters:
  • etarho_edge – face-centered \(\eta_\rho\)

  • etarho_cell – cell-centered \(\eta_\rho\)

  • etarho_flux\(\eta_\rho\) flux

void MakeEtarhoSphr(const amrex::Vector<amrex::MultiFab> &scal_old, const amrex::Vector<amrex::MultiFab> &scal_new, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac)
void MakeEtarhoPlanar(const amrex::Vector<amrex::MultiFab> &scal_old, const amrex::Vector<amrex::MultiFab> &scal_new, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac)
void MakeGravCell(BaseState<amrex::Real> &grav_cell, const BaseState<amrex::Real> &rho0_s)
void MakeGravEdge(BaseState<amrex::Real> &grav_edge_state, const BaseState<amrex::Real> &rho0_state)
void MakePsiPlanar()
void MakePsiSphr(const BaseState<amrex::Real> &gamma1bar, const BaseState<amrex::Real> &p0_avg, const BaseState<amrex::Real> &Sbar_in)
void MakePsiIrreg(const BaseState<amrex::Real> &grav_cell)
void Make_S_cc(amrex::Vector<amrex::MultiFab> &S_cc, amrex::Vector<amrex::MultiFab> &delta_gamma1_term, amrex::Vector<amrex::MultiFab> &delta_gamma1, const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &u, const amrex::Vector<amrex::MultiFab> &rho_omegadot, const amrex::Vector<amrex::MultiFab> &rho_Hnuc, const amrex::Vector<amrex::MultiFab> &rho_Hext, const amrex::Vector<amrex::MultiFab> &thermal, const BaseState<amrex::Real> &p0_s, const BaseState<amrex::Real> &gamma1bar, BaseState<amrex::Real> &delta_gamma1_termbar)

Compute S at cell-centers.

void MakeRHCCforNodalProj(amrex::Vector<amrex::MultiFab> &rhcc, const amrex::Vector<amrex::MultiFab> &S_cc, const BaseState<amrex::Real> &Sbar, const BaseState<amrex::Real> &beta0, const amrex::Vector<amrex::MultiFab> &delta_gamma1_term)

Compute rhcc = beta0*(S_cc-Sbar) + beta0*delta_chi

void MakeRHCCforMacProj(amrex::Vector<amrex::MultiFab> &rhcc, const BaseState<amrex::Real> &rho0, const amrex::Vector<amrex::MultiFab> &S_cc, const BaseState<amrex::Real> &Sbar, const BaseState<amrex::Real> &beta0, const amrex::Vector<amrex::MultiFab> &delta_gamma1_term, const BaseState<amrex::Real> &gamma1bar, const BaseState<amrex::Real> &p0, const amrex::Vector<amrex::MultiFab> &delta_p_term, amrex::Vector<amrex::MultiFab> &delta_chi, const bool is_predictor)

Compute rhcc = beta0*(S_cc-Sbar) + beta0*delta_chi

void CorrectRHCCforNodalProj(amrex::Vector<amrex::MultiFab> &rhcc, const BaseState<amrex::Real> &rho0, const BaseState<amrex::Real> &beta0, const BaseState<amrex::Real> &gamma1bar, const BaseState<amrex::Real> &p0, const amrex::Vector<amrex::MultiFab> &delta_p_term)
void MakeUtrans(const amrex::Vector<amrex::MultiFab> &utilde, const amrex::Vector<amrex::MultiFab> &ufull, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &utrans, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac)

Create utrans, the transverse velocity

Parameters:
  • utilde – perturbed velocity

  • ufull – full velocity

  • utrans – transverse velocity

  • w0mac – MAC base-state velocity

void Makew0(const BaseState<amrex::Real> &w0_old, BaseState<amrex::Real> &w0_force, const BaseState<amrex::Real> &Sbar_in, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &p0_old_in, const BaseState<amrex::Real> &p0_new_in, const BaseState<amrex::Real> &gamma1bar_old_in, const BaseState<amrex::Real> &gamma1bar_new_in, const BaseState<amrex::Real> &p0_minus_peosbar, const amrex::Real dt_in, const amrex::Real dtold_in, const bool is_predictor)
void Makew0PlanarVarg(const BaseState<amrex::Real> &w0_old, BaseState<amrex::Real> &w0_force, const BaseState<amrex::Real> &Sbar_in, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &p0_old_in, const BaseState<amrex::Real> &p0_new_in, const BaseState<amrex::Real> &gamma1bar_old_in, const BaseState<amrex::Real> &gamma1bar_new_in, const BaseState<amrex::Real> &p0_minus_peosbar, const amrex::Real dt_in, const amrex::Real dtold_in)
void Makew0Planar(const BaseState<amrex::Real> &w0_old, BaseState<amrex::Real> &w0_force, const BaseState<amrex::Real> &Sbar_in, const BaseState<amrex::Real> &p0_old_in, const BaseState<amrex::Real> &p0_new_in, const BaseState<amrex::Real> &gamma1bar_old_in, const BaseState<amrex::Real> &gamma1bar_new_in, const BaseState<amrex::Real> &p0_minus_peosbar, const amrex::Real dt_in, const amrex::Real dtold_in, const bool is_predictor)
void Makew0SphrIrreg(const BaseState<amrex::Real> &w0_old, BaseState<amrex::Real> &w0_force, const BaseState<amrex::Real> &Sbar_in, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &p0_old_in, const BaseState<amrex::Real> &p0_new_in, const BaseState<amrex::Real> &gamma1bar_old_in, const BaseState<amrex::Real> &gamma1bar_new_in, const BaseState<amrex::Real> &p0_minus_peosbar, const amrex::Real dt_in, const amrex::Real dtold_in)
void Makew0Sphr(const BaseState<amrex::Real> &w0_old, BaseState<amrex::Real> &w0_force, const BaseState<amrex::Real> &Sbar_in, const BaseState<amrex::Real> &rho0_old_in, const BaseState<amrex::Real> &rho0_new_in, const BaseState<amrex::Real> &p0_old_in, const BaseState<amrex::Real> &p0_new_in, const BaseState<amrex::Real> &gamma1bar_old_in, const BaseState<amrex::Real> &gamma1bar_new_in, const BaseState<amrex::Real> &p0_minus_peosbar, const amrex::Real dt_in, const amrex::Real dtold_in)
void Tridiag(const BaseStateArray<Real> &a, const BaseStateArray<Real> &b, const BaseStateArray<Real> &c, const BaseStateArray<Real> &r, const BaseStateArray<Real> &u, const int n)
void ProlongBasetoUniform(const BaseState<amrex::Real> &base_ml_s, BaseState<amrex::Real> &base_fine_s) const
void NodalProj(int proj_type, amrex::Vector<amrex::MultiFab> &rhcc, int istep_divu_iter = 0, bool sdc_off = true)

Perform a nodal projection.

void CreateUvecForProj(int proj_type, amrex::Vector<amrex::MultiFab> &Vproj, const amrex::Vector<amrex::MultiFab> &sig)

Fill in Vproj

initial_projection_comp: Utilde^0 -- uold divu_iters_comp: Utilde^0 -- uold pressure_iters_comp: (Utilde^n+1,* - Utilde^n)/dt -- (unew-uold)/dt regular_timestep_comp: (Utilde^n+1,* + dt*gpi/rhohalf) -- unew + dt*gpi/rhohalf

Parameters:

sig – contains rhohalf if proj_type == regular_timestep_comp

void SetBoundaryVelocity(amrex::Vector<amrex::MultiFab> &vel)
void ComputeGradPhi(amrex::Vector<amrex::MultiFab> &phi, amrex::Vector<amrex::MultiFab> &gphi)

Given a nodal phi, compute \(\nabla(\phi)\) at cell centers.

void MakePiCC(const amrex::Vector<amrex::MultiFab> &beta0_cart)

Average nodal pi to cell-centers and put in the Pi component of snew

void PlotFileName(const int lev, std::string *plotfilename)

Get plotfile name.

amrex::Vector<const amrex::MultiFab*> PlotFileMF(const int nPlot, const amrex::Real t_in, const amrex::Real dt_in, const amrex::Vector<amrex::MultiFab> &rho0_cart, const amrex::Vector<amrex::MultiFab> &rhoh0_cart, const amrex::Vector<amrex::MultiFab> &p0_cart, const amrex::Vector<amrex::MultiFab> &gamma1bar_cart, const amrex::Vector<amrex::MultiFab> &u_in, amrex::Vector<amrex::MultiFab> &s_in, const BaseState<amrex::Real> &p0_in, const BaseState<amrex::Real> &gamma1bar_in, const amrex::Vector<amrex::MultiFab> &S_cc_in)

Put together an array of multifabs for writing.

amrex::Vector<const amrex::MultiFab*> SmallPlotFileMF(const int nPlot, const int nSmallPlot, amrex::Vector<const amrex::MultiFab*> mf, const amrex::Vector<std::string> &varnames, const amrex::Vector<std::string> &small_plot_varnames)
amrex::Vector<std::string> PlotFileVarNames(int *nPlot) const

Set plotfile variables names.

amrex::Vector<std::string> SmallPlotFileVarNames(int *nPlot, amrex::Vector<std::string> varnames) const

Set small plotfile variables names.

void WriteSmallPlotFile(const int step, const amrex::Real t_in, const amrex::Real dt_in, const BaseState<amrex::Real> &rho0_in, const BaseState<amrex::Real> &rhoh0_in, const BaseState<amrex::Real> &p0_in, const BaseState<amrex::Real> &gamma1bar_in, const amrex::Vector<amrex::MultiFab> &u_in, amrex::Vector<amrex::MultiFab> &s_in, const amrex::Vector<amrex::MultiFab> &S_cc_in)

Write a small plotfile to disk.

void WritePlotFile(const int step, const amrex::Real t_in, const amrex::Real dt_in, const BaseState<amrex::Real> &rho0_in, const BaseState<amrex::Real> &rhoh0_in, const BaseState<amrex::Real> &p0_in, const BaseState<amrex::Real> &gamma1bar_in, const amrex::Vector<amrex::MultiFab> &u_in, amrex::Vector<amrex::MultiFab> &s_in, const amrex::Vector<amrex::MultiFab> &S_cc_in, const bool is_small = false)

Write plotfile to disk.

void WriteJobInfo(const std::string &dir) const
void MakeMagvel(const amrex::Vector<amrex::MultiFab> &vel, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, amrex::Vector<amrex::MultiFab> &magvel)

Calculate the magnitude of the velocity.

void MakeVelrc(const amrex::Vector<amrex::MultiFab> &vel, const amrex::Vector<amrex::MultiFab> &w0rcart, amrex::Vector<amrex::MultiFab> &rad_vel, amrex::Vector<amrex::MultiFab> &circ_vel)

Calculate the radial and circular components of the velocity.

void MakeAdExcess(const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<amrex::MultiFab> &ad_excess)

Calculate the adiabatic excess.

void MakeGrav(const BaseState<amrex::Real> &rho0, amrex::Vector<amrex::MultiFab> &grav)

Calculate the gravitational acceleration.

void MakeVorticity(const amrex::Vector<amrex::MultiFab> &vel, amrex::Vector<amrex::MultiFab> &vorticity)

Calculate the vorticity.

void MakeDeltaGamma(const amrex::Vector<amrex::MultiFab> &state, const BaseState<amrex::Real> &p0, const amrex::Vector<amrex::MultiFab> &p0_cart, const BaseState<amrex::Real> &gamma1bar, const amrex::Vector<amrex::MultiFab> &gamma1bar_cart, amrex::Vector<amrex::MultiFab> &deltagamma)

Calculate deltagamma

void MakeEntropy(const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<amrex::MultiFab> &entropy)

Calculate the entropy.

void MakeDivw0(const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, amrex::Vector<amrex::MultiFab> &divw0)

Calculate the divergenve of the base state velocity.

void MakePiDivu(const amrex::Vector<amrex::MultiFab> &vel, const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<amrex::MultiFab> &pidivu)

Calculate pi times the divergence of the velocity.

void MakeAbar(const amrex::Vector<amrex::MultiFab> &state, amrex::Vector<amrex::MultiFab> &abar)

Mass fractions of the species.

void PPM(const amrex::Box &bx, amrex::Array4<const amrex::Real> const s, amrex::Array4<const amrex::Real> const u, amrex::Array4<const amrex::Real> const v, amrex::Array4<const amrex::Real> const w, amrex::Array4<amrex::Real> const Ip, amrex::Array4<amrex::Real> const Im, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx, const bool is_umac, const int comp, const int bccomp) const
void React(const amrex::Vector<amrex::MultiFab> &s_in, amrex::Vector<amrex::MultiFab> &s_out, amrex::Vector<amrex::MultiFab> &rho_Hext, amrex::Vector<amrex::MultiFab> &rho_omegadot, amrex::Vector<amrex::MultiFab> &rho_Hnuc, const BaseState<amrex::Real> &p0, const amrex::Real dt_in, const amrex::Real time_in)

Compute heating term, rho_Hext, then react the state over dt_react and update rho_omegadot, rho_Hnuc

void Burner(const amrex::Vector<amrex::MultiFab> &s_in, amrex::Vector<amrex::MultiFab> &s_out, const amrex::Vector<amrex::MultiFab> &rho_Hext, amrex::Vector<amrex::MultiFab> &rho_omegadot, amrex::Vector<amrex::MultiFab> &rho_Hnuc, const BaseState<amrex::Real> &p0, const amrex::Real dt_in, const amrex::Real time_in)
void MakeIntraCoeffs(const amrex::Vector<amrex::MultiFab> &scal1, const amrex::Vector<amrex::MultiFab> &scal2, amrex::Vector<amrex::MultiFab> &cp, amrex::Vector<amrex::MultiFab> &xi)
void MakeHeating(amrex::Vector<amrex::MultiFab> &rho_Hext, const amrex::Vector<amrex::MultiFab> &scal)
void Regrid()

Check to see if we need to regrid, then regrid.

void TagArray()

Set tagging array to include buffer zones for multilevel.

void RegridBaseState(BaseState<amrex::Real> &base_s, const bool is_edge = false)

Regrid base state variables (ex. psi, etarho, rho0, etc.)

We copy the coarsest level only, interpolate to all the other levels and then copy the valid data from the old arrays onto the new.

void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ng) override

Tag all cells for refinement

Overrides the pure virtual function in AmrCore

void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override

Within a call to AmrCore::regrid, this function fills in data at a level that existed before, using pre-existing fine and interpolated coarse data

overrides the pure virtual function in AmrCore

void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override

Within a call to AmrCore::regrid, this function fills in data at a level that did NOT exist before, using interpolated coarse data

overrides the pure virtual function in AmrCore

void ClearLevel(int lev) override

Within a call to AmrCore::regrid, this function deletes all data at a level of refinement that is no longer needed

overrides the pure virtual function in AmrCore

void TfromRhoH(amrex::Vector<amrex::MultiFab> &scal, const BaseState<amrex::Real> &p0)

Calculate the temperature given the density and the enthalpy

Parameters:
  • scal – scalars

  • p0 – base state pressure

void TfromRhoP(amrex::Vector<amrex::MultiFab> &scal, const BaseState<amrex::Real> &p0, const bool updateRhoH = false)

Calculate the temperature given the density and the pressure

Parameters:
  • scal – scalars

  • p0 – base state pressure

void PfromRhoH(const amrex::Vector<amrex::MultiFab> &state, const amrex::Vector<amrex::MultiFab> &s_old, amrex::Vector<amrex::MultiFab> &peos)

Calculate the pressure given the density and the enthalpy

Parameters:
  • state – scalars

  • s_old – scalars at old time step

  • peos – pressure calculated from the equation of state

void MachfromRhoH(const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &vel, const BaseState<amrex::Real> &p0, const amrex::Vector<amrex::MultiFab> &w0cart, amrex::Vector<amrex::MultiFab> &mach)

Calculate the Mach number given the density and the enthalpy

Parameters:
  • scal – scalars

  • vel – velocity

  • p0 – base state pressure

  • mach – Mach number

void CsfromRhoH(const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &p0_cart, amrex::Vector<amrex::MultiFab> &cs)

Calculate the sound speed given the density and the enthalpy

Parameters:
  • scal – scalars

  • p0 – base state pressure

  • p0cart – base state pressure on cartesian grid

  • cs – sound speed

void HfromRhoTedge(amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &sedge, const BaseState<amrex::Real> &rho0_edge_old, const BaseState<amrex::Real> &rhoh0_edge_old, const BaseState<amrex::Real> &rho0_edge_new, const BaseState<amrex::Real> &rhoh0_edge_new)
void RotationInit()
void ExternInit()
void ReadParameters()

Read in C++ parameters from inputs file.

void VariableSetup()

Define variable mappings (Rho, RhoH, …, Nscal, etc.)

void BCSetup()

Set up BCRec definitions for BC types.

void Slopex(const amrex::Box &bx, amrex::Array4<amrex::Real> const s, amrex::Array4<amrex::Real> const slx, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, int ncomp, int bc_start_comp)
void Slopey(const amrex::Box &bx, amrex::Array4<amrex::Real> const s, amrex::Array4<amrex::Real> const sly, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, int ncomp, int bc_start_comp)
void Slopez(const amrex::Box &bx, amrex::Array4<amrex::Real> const s, amrex::Array4<amrex::Real> const slz, const amrex::Box &domainBox, const amrex::Vector<amrex::BCRec> &bcs, int ncomp, int bc_start_comp)
void SpongeInit(const BaseState<amrex::Real> &rho0_s)
void MakeSponge(amrex::Vector<amrex::MultiFab> &sponge)
void RetagArray(const amrex::Box &bx, const int lev)
void TagBoxes(amrex::TagBoxArray &tags, const amrex::MFIter &mfi, const int lev, const amrex::Real time)
void StateError(amrex::TagBoxArray &tags, const amrex::MultiFab &state_mf, const amrex::MFIter &mfi, const int lev, const amrex::Real time)
void MakeExplicitThermal(amrex::Vector<amrex::MultiFab> &thermal, const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &Tcoeff, const amrex::Vector<amrex::MultiFab> &hcoeff, const amrex::Vector<amrex::MultiFab> &Xkcoeff, const amrex::Vector<amrex::MultiFab> &pcoeff, const BaseState<amrex::Real> &p0, int temp_formulation)

Compute the quantity: thermal = del dot kappa grad T

  • if temp_diffusion_formulation = 1, then we compute this directly.

  • if temp_diffusion_formulation = 2, then we compute the algebraically equivalent form with grad h - grad X_k - grad p_0 formulation

Parameters:
  • thermal – thermal term

  • scal – scalars

  • Tcoeff – temperature coefficient

  • hcoeff – enthalpy coefficient

  • Xkcoeff – species coefficients

  • pcoeff – pressure coefficient

  • p0 – base state pressure

void ApplyThermal(amrex::MLABecLaplacian &mlabec, amrex::Vector<amrex::MultiFab> &thermalout, const amrex::Vector<amrex::MultiFab> &coeff, amrex::Vector<amrex::MultiFab> &phi, const amrex::Vector<amrex::BCRec> &bcs, int bccomp)

Use apply() to construct the form of the conduction term. apply() forms the generic quantity:

(alpha * A - beta * div B grad) phi = RHS

void MakeThermalCoeffs(const amrex::Vector<amrex::MultiFab> &scal, amrex::Vector<amrex::MultiFab> &Tcoeff, amrex::Vector<amrex::MultiFab> &hcoeff, amrex::Vector<amrex::MultiFab> &Xkcoeff, amrex::Vector<amrex::MultiFab> &pcoeff)

create the coefficients for grad{T}, grad{h}, grad{X_k}, and grad{p_0} for the thermal diffusion term in the enthalpy equation.

note: we explicitly fill the ghostcells by looping over them directly

Parameters:
  • scal – scalars

  • Tcoeff – temperature coefficient

  • hcoeff – enthalpy coefficient

  • Xkcoeff – species coefficients

  • pcoeff – pressure coefficient

void ThermalConduct(const amrex::Vector<amrex::MultiFab> &s1, amrex::Vector<amrex::MultiFab> &s2, const amrex::Vector<amrex::MultiFab> &hcoeff1, const amrex::Vector<amrex::MultiFab> &Xkcoeff1, const amrex::Vector<amrex::MultiFab> &pcoeff1, const amrex::Vector<amrex::MultiFab> &hcoeff2, const amrex::Vector<amrex::MultiFab> &Xkcoeff2, const amrex::Vector<amrex::MultiFab> &pcoeff2)

ThermalConduct implements thermal diffusion in the enthalpy equation. This is an implicit solve, using the multigrid solver. This updates the enthalpy only.

void MakeExplicitThermalHterm(amrex::Vector<amrex::MultiFab> &thermal, const amrex::Vector<amrex::MultiFab> &scal, const amrex::Vector<amrex::MultiFab> &hcoeff)
void VelocityAdvance(const amrex::Vector<amrex::MultiFab> &rhohalf, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const amrex::Vector<amrex::MultiFab> &w0_force_cart, const BaseState<amrex::Real> &rho0_nph, const BaseState<amrex::Real> &grav_cell_nph, const amrex::Vector<amrex::MultiFab> &sponge)

Advance the velocity.

void VelPred(amrex::Vector<amrex::MultiFab> &utilde, const amrex::Vector<amrex::MultiFab> &ufull, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &utrans, amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &umac, const amrex::Vector<std::array<amrex::MultiFab, AMREX_SPACEDIM>> &w0mac, const amrex::Vector<amrex::MultiFab> &force)

Predict the normal velocities to the interfaces. We don’t care about the transverse velocities here. The prediction is done piecewise linear (for now)

Parameters:
  • utilde – perturbed velocity

  • ufull – full velocity

  • utrans – transverse velocity

  • umac – MAC velocity

  • w0mac – MAC base-state velocity

  • force – velocity force

void VelPredInterface(const amrex::MFIter &mfi, amrex::Array4<const amrex::Real> const utilde, amrex::Array4<const amrex::Real> const ufull, amrex::Array4<const amrex::Real> const utrans, amrex::Array4<const amrex::Real> const vtrans, amrex::Array4<const amrex::Real> const Imu, amrex::Array4<const amrex::Real> const Ipu, amrex::Array4<const amrex::Real> const Imv, amrex::Array4<const amrex::Real> const Ipv, amrex::Array4<amrex::Real> const ulx, amrex::Array4<amrex::Real> const urx, amrex::Array4<amrex::Real> const uimhx, amrex::Array4<amrex::Real> const uly, amrex::Array4<amrex::Real> const ury, amrex::Array4<amrex::Real> const uimhy, const amrex::Box &domainBox, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx)
void VelPredVelocities(const amrex::MFIter &mfi, amrex::Array4<const amrex::Real> const utilde, amrex::Array4<const amrex::Real> const utrans, amrex::Array4<const amrex::Real> const vtrans, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<const amrex::Real> const Imfx, amrex::Array4<const amrex::Real> const Ipfx, amrex::Array4<const amrex::Real> const Imfy, amrex::Array4<const amrex::Real> const Ipfy, amrex::Array4<const amrex::Real> const ulx, amrex::Array4<const amrex::Real> const urx, amrex::Array4<const amrex::Real> const uimhx, amrex::Array4<const amrex::Real> const uly, amrex::Array4<const amrex::Real> const ury, amrex::Array4<const amrex::Real> const uimhy, amrex::Array4<const amrex::Real> const force, amrex::Array4<const amrex::Real> const w0_cart_in, const amrex::Box &domainBox, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx)
void VelPredInterface(const amrex::MFIter &mfi, amrex::Array4<const amrex::Real> const utilde, amrex::Array4<const amrex::Real> const ufull, amrex::Array4<const amrex::Real> const utrans, amrex::Array4<const amrex::Real> const vtrans, amrex::Array4<const amrex::Real> const wtrans, amrex::Array4<const amrex::Real> const Imu, amrex::Array4<const amrex::Real> const Ipu, amrex::Array4<const amrex::Real> const Imv, amrex::Array4<const amrex::Real> const Ipv, amrex::Array4<const amrex::Real> const Imw, amrex::Array4<const amrex::Real> const Ipw, amrex::Array4<amrex::Real> const ulx, amrex::Array4<amrex::Real> const urx, amrex::Array4<amrex::Real> const uimhx, amrex::Array4<amrex::Real> const uly, amrex::Array4<amrex::Real> const ury, amrex::Array4<amrex::Real> const uimhy, amrex::Array4<amrex::Real> const ulz, amrex::Array4<amrex::Real> const urz, amrex::Array4<amrex::Real> const uimhz, const amrex::Box &domainBox, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx)
void VelPredTransverse(const amrex::MFIter &mfi, amrex::Array4<const amrex::Real> const utilde, amrex::Array4<const amrex::Real> const utrans, amrex::Array4<const amrex::Real> const vtrans, amrex::Array4<const amrex::Real> const wtrans, amrex::Array4<const amrex::Real> const ulx, amrex::Array4<const amrex::Real> const urx, amrex::Array4<const amrex::Real> const uimhx, amrex::Array4<const amrex::Real> const uly, amrex::Array4<const amrex::Real> const ury, amrex::Array4<const amrex::Real> const uimhy, amrex::Array4<const amrex::Real> const ulz, amrex::Array4<const amrex::Real> const urz, amrex::Array4<const amrex::Real> const uimhz, amrex::Array4<amrex::Real> const uimhyz, amrex::Array4<amrex::Real> const uimhzy, amrex::Array4<amrex::Real> const vimhxz, amrex::Array4<amrex::Real> const vimhzx, amrex::Array4<amrex::Real> const wimhxy, amrex::Array4<amrex::Real> const wimhyx, const amrex::Box &domainBox, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx)
void VelPredVelocities(const amrex::MFIter &mfi, amrex::Array4<const amrex::Real> const utilde, amrex::Array4<const amrex::Real> const utrans, amrex::Array4<const amrex::Real> const vtrans, amrex::Array4<const amrex::Real> const wtrans, amrex::Array4<amrex::Real> const umac, amrex::Array4<amrex::Real> const vmac, amrex::Array4<amrex::Real> const wmac, amrex::Array4<const amrex::Real> const w0macx, amrex::Array4<const amrex::Real> const w0macy, amrex::Array4<const amrex::Real> const w0macz, amrex::Array4<const amrex::Real> const Imfx, amrex::Array4<const amrex::Real> const Ipfx, amrex::Array4<const amrex::Real> const Imfy, amrex::Array4<const amrex::Real> const Ipfy, amrex::Array4<const amrex::Real> const Imfz, amrex::Array4<const amrex::Real> const Ipfz, amrex::Array4<const amrex::Real> const ulx, amrex::Array4<const amrex::Real> const urx, amrex::Array4<const amrex::Real> const uly, amrex::Array4<const amrex::Real> const ury, amrex::Array4<const amrex::Real> const ulz, amrex::Array4<const amrex::Real> const urz, amrex::Array4<const amrex::Real> const uimhyz, amrex::Array4<const amrex::Real> const uimhzy, amrex::Array4<const amrex::Real> const vimhxz, amrex::Array4<const amrex::Real> const vimhzx, amrex::Array4<const amrex::Real> const wimhxy, amrex::Array4<const amrex::Real> const wimhyx, amrex::Array4<const amrex::Real> const force, amrex::Array4<const amrex::Real> const w0_cart_in, const amrex::Box &domainBox, const amrex::GpuArray<Real, AMREX_SPACEDIM> dx)

Public Members

int istep

which step?

int start_step
amrex::Real t_new
amrex::Real t_old
amrex::Real dt
amrex::Real dtold
int ng_adv

number of ghost cells needed for hyperbolic step

ModelParser input_model
SimpleLog log_file
amrex::Vector<amrex::MultiFab> sold
amrex::Vector<amrex::MultiFab> snew
amrex::Vector<amrex::MultiFab> uold
amrex::Vector<amrex::MultiFab> unew
amrex::Vector<amrex::MultiFab> S_cc_old
amrex::Vector<amrex::MultiFab> S_cc_new
amrex::Vector<amrex::MultiFab> gpi
amrex::Vector<amrex::MultiFab> dSdt
amrex::Vector<amrex::MultiFab> pi
amrex::Vector<amrex::MultiFab> w0_cart

this doesn’t have to be persistent, but we make it so that we avoid continually creating and filling temporaries saves on some flops and data movement (GPU)

amrex::Vector<amrex::MultiFab> rhcc_for_nodalproj

this only needs to persist leading into the initial pressure iters since we project (beta0^nph S^1 - beta0 S^0) / dt during a regular time step we overwrite this

amrex::Vector<amrex::MultiFab> normal

spherical only - we make this persistent in that we only have to rebuild and fill this after regridding

amrex::Vector<amrex::iMultiFab> cell_cc_to_r
IntVector phys_bc

stores domain boundary conditions. These muse be vectors (rather than arrays) so we can ParmParse them

amrex::Vector<amrex::BCRec> bcs_s

Boundary condition objects needed for FillPatch routines. This is essentially an array (over components) of 2*DIM integer arrays storing the physical boundary condition types at the lo/hi walls in each direction

amrex::Vector<amrex::BCRec> bcs_u
amrex::Vector<amrex::BCRec> bcs_f
BaseState<amrex::Real> s0_init
BaseState<amrex::Real> p0_init
BaseState<amrex::Real> rho0_old
BaseState<amrex::Real> rho0_new
BaseState<amrex::Real> rhoh0_old
BaseState<amrex::Real> rhoh0_new
BaseState<amrex::Real> p0_old
BaseState<amrex::Real> p0_new
BaseState<amrex::Real> tempbar
BaseState<amrex::Real> tempbar_init
BaseState<amrex::Real> beta0_old
BaseState<amrex::Real> beta0_new
BaseState<amrex::Real> gamma1bar_old
BaseState<amrex::Real> gamma1bar_new
BaseState<amrex::Real> etarho_cc
BaseState<amrex::Real> psi
BaseState<amrex::Real> grav_cell_old
BaseState<amrex::Real> grav_cell_new
BaseState<amrex::Real> p0_nm1

p0^{n-1} needed to compute d(p0)/dt for nonuniform grid spacing

BaseState<amrex::Real> beta0_nm1
BaseState<amrex::Real> w0
BaseState<amrex::Real> etarho_ec
amrex::Real rel_eps
IntVector tag_array

array of tagged boxes (planar)

BaseStateGeometry base_geom

contains base state geometry variables

amrex::Vector<amrex::Real> diagfile1_data
amrex::Vector<amrex::Real> diagfile2_data
amrex::Vector<amrex::Real> diagfile3_data
amrex::GpuArray<Real, 3> center
Real sin_theta
Real cos_theta
Real omega
Real r_sp
Real r_md
Real r_tp
Real r_sp_outer
Real r_tp_outer
Real sponge_start_density

the sponge_start_density should be the density below which the sponge first turns on. Different problems may compute this in different ways (i.e. not using sponge_center_density and sponge_start_factor), so we provide this public module variable to ensure that the rest of the code always knows at what density the sponge begins.

amrex::Vector<std::unique_ptr<amrex::FluxRegister>> flux_reg_s

Stores fluxes at coarse-fine interface for synchronization. This will be sized max_level+1

NOTE: the flux register associated with flux_reg[lev] is associated with the lev/lev-1 interface (and has grid spacing associated with lev-1) therefore flux_reg[0] and flux_reg[max_level] are never actually used in the reflux operation

Real p0bdot

Rate of change of pressure at the lower boundary. Used to modify w0 and p0 in the case of a closed upper boundary (planar).

Real p0b

Public Static Functions

static amrex::Real getCPUTime()

Get the amount of CPU time used. This will persist after restarts.

static void WriteBuildInfo()

Dump build info.

static void ScalarFill(const amrex::Array4<amrex::Real> &scal, const amrex::Box &bx, const amrex::Box &domainBox, const amrex::Real *dx, const amrex::BCRec bcs, const amrex::Real *gridlo, const int comp)
static void VelFill(const amrex::Array4<amrex::Real> &vel, const amrex::Box &bx, const amrex::Box &domainBox, const amrex::Real *dx, const amrex::BCRec bcs, const amrex::Real *gridlo, const int comp = 0)
static void FillExtBC(const amrex::Array4<amrex::Real> &q, const amrex::Box &bx, const amrex::Box &domainBox, const amrex::Real *dx, const amrex::BCRec bcs, const int comp, const bool is_vel = false)

Public Static Attributes

static amrex::Real previousCPUTimeUsed
static amrex::Real startCPUTime
static int ng_s
static amrex::IntVect nodal_flag
static amrex::IntVect nodal_flag_x
static amrex::IntVect nodal_flag_y
static amrex::IntVect nodal_flag_z
static constexpr amrex::Real Gconst = 6.67428e-8