Class Castro
-
class Castro : public amrex::AmrLevel
AmrLevel
-derived class for hyperbolic conservation equations for stellar mediaPublic Functions
-
Castro()
Default constructor. Builds invalid object.
-
Castro(amrex::Amr &papa, int lev, const amrex::Geometry &level_geom, const amrex::BoxArray &bl, const amrex::DistributionMapping &dm, amrex::Real time)
The basic constructor.
- Parameters:
papa – parent
amrex::Amr
objectlev – which level are we on
level_geom – level geometry
bl –
BoxArray
objectdm – the mapping
time – current time
-
~Castro() override
The destructor.
-
void restart(amrex::Amr &papa, istream &is, bool bReadSpecial = false) override
Restart from a checkpoint file.
- Parameters:
papa – parent
amrex::Amr
objectis –
bReadSpecial –
-
void set_state_in_checkpoint(amrex::Vector<int> &state_in_checkpoint) override
This is called only when we restart from an old checkpoint. Sets the Vector
state_in_checkpoint
depending on the input version and type.- Parameters:
state_in_checkpoint – Vector of integer flags describing the state
-
void checkPoint(const std::string &dir, std::ostream &os, amrex::VisMF::How how, bool dump_old) override
Call
amrex::AmrLevel::checkPoint
and then add radiation info- Parameters:
dir – Directory to store checkpoint in
os –
std::ostream
objecthow –
V
isMF::How`` objectdump_old – boolean
-
std::string thePlotFileType() const override
A string written as the first item in writePlotFile() at level zero. It is so we can distinguish between different types of plot files. For Castro it has the form: Castro-Vnnn
-
void setPlotVariables() override
-
void writePlotFile(const std::string &dir, ostream &os, amrex::VisMF::How how) override
Write a plotfile to specified directory.
- Parameters:
dir – directory to write to
os – output stream
how –
VisMF::How
object
-
void writeSmallPlotFile(const std::string &dir, ostream &os, amrex::VisMF::How how) override
- Parameters:
dir – directory to write to
os – output stream
how –
VisMF::How
object
-
void plotFileOutput(const std::string &dir, ostream &os, amrex::VisMF::How how, const int is_small)
- Parameters:
dir – directory to write to
os – output stream
how –
VisMF::How
objectis_small – integer
-
void writeJobInfo(const std::string &dir, const amrex::Real io_time)
Write job info to file
- Parameters:
dir – directory to write to
io_time – how long it took to generate the plotfile/checkpoint
-
void initData() override
Initialize grid data at problem start-up.
-
void init_particles()
Initialize particle locations and velocities (and strengths if relevant)
-
void ParticleCheckPoint(const std::string &dir)
Write particles in checkpoint directories
- Parameters:
dir – directory containing checkpoint
-
void ParticlePlotFile(const std::string &dir)
Write particles in plotfile directories
- Parameters:
dir – directory containing plotfile
-
void ParticlePostRestart(const std::string &restart_file)
How to initialize at restart
- Parameters:
restart_file – path to restart file
-
std::unique_ptr<amrex::MultiFab> ParticleDerive(const std::string &name, amrex::Real time, int ngrow)
Derived quantities associated with particles
- Parameters:
name – name of particle variable to derive
time – current time
ngrow – Number of ghost cells
-
void TimestampParticles(int ngrow)
Timestamp particles
- Parameters:
ngrow –
-
void advance_particles(int iteration, amrex::Real time, amrex::Real dt)
Advance the particles by dt
- Parameters:
iteration – where we are in the current AMR subcycle
time – current time
dt – timestep
-
void MAESTRO_init()
-
amrex::MultiFab &Area(int dir)
Area of the multifab with surface normal along
dir
- Parameters:
dir – direction of normal
-
amrex::MultiFab &Fluxes(int dir)
Return the n’th fluxes MultiFab.
- Parameters:
dir – direction in which to find the fluxes along
-
void setTimeLevel(amrex::Real time, amrex::Real dt_old, amrex::Real dt_new) override
Set time levels of state data.
- Parameters:
time – current time
dt_old – old timestep
dt_new – new timestep
-
void init(amrex::AmrLevel &old) override
Initialize data on this level from another Castro (during regrid).
- Parameters:
old – amrex::AmrLevel object to use to initialize data
-
void init() override
Initialize data on this level after regridding if old level did not previously exist
-
int okToContinue() override
Proceed with next timestep?
-
amrex::Real advance(amrex::Real time, amrex::Real dt, int iteration, int ncycle) override
Advance grids at this level in time.
The main driver for a single level. This will do either the SDC algorithm or the Strang-split reactions algorithm.
- Parameters:
time – the current simulation time
dt – the timestep to advance (e.g., go from time to time + dt)
iteration – where we are in the current AMR subcycle. Each level will take a number of steps to reach the final time of the coarser level below it. This counter starts at 1
ncycle – the number of subcycles at this level
-
advance_status do_advance_ctu(amrex::Real time, amrex::Real dt)
This routine will advance the old state data (called
S_old
here) to the new time, for a single level. The new data is calledS_new
here. The update includes reactions (if we are not doing SDC), hydro, and the source terms. The advance is doing using the CTU (corner transport upwind) approach. We also use this interface for the simplified SDC update.- Parameters:
time – the current simulation time
dt – the timestep to advance (e.g., go from time to time + dt)
-
void save_data_for_retry()
Save a copy of the old state data in case for the purposes of a retry.
-
bool retry_advance_ctu(amrex::Real dt, const advance_status &status)
Should we retry advancing the simulation? By default, we don’t do a retry unless the criteria are violated. This is used only for the CTU advance.
- Parameters:
dt – the timestep to advance (e.g., go from time to time + dt)
advance_status – was the advance successful?
-
amrex::Real subcycle_advance_ctu(amrex::Real time, amrex::Real dt, int amr_iteration, int amr_ncycle)
Subcyles until we’ve reached the target time,
time
+dt
. The last timestep will be shortened if needed so that it does not overshoot the ending time. This is used only for the CTU advance.- Parameters:
time – the current simulation time
dt – the timestep to advance (e.g., go from time to time + dt)
amr_iteration – where we are in the current AMR subcycle. Each level will take a number of steps to reach the final time of the coarser level below it. This counter starts at 1
amr_ncycle – the number of subcycles at this level
-
void initialize_advance(amrex::Real time, amrex::Real dt, int amr_iteration)
This is run at the start of an ::advance on a single level. It sets up all the necessary arrays and data containers, ensures the data is valid and makes copies of the MultiFabs in the old and new state in case a retry is performed.
- Parameters:
time – the current simulation time
dt – the timestep to advance (e.g., go from time to time + dt)
amr_iteration – where we are in the current AMR subcycle. Each level will take a number of steps to reach the final time of the coarser level below it. This counter starts at 1
-
void finalize_advance()
If required, adds material lost in timestep to cumulative losses, does refluxing, and clears arrays/data containers initialized at start of the level’s advance.
-
advance_status initialize_do_advance(amrex::Real time, amrex::Real dt)
Performed at the start of ::do_advance, resets flags and ensures required ghost zones are filled.
- Parameters:
time – the current simulation time
dt – the timestep to advance
-
advance_status finalize_do_advance(amrex::Real time, amrex::Real dt)
Performed at the end of ::do_advance, cleans up temporary data created by ::initialize_do_advance
- Parameters:
time – the current simulation time
dt – the timestep to advance
-
void computeInitialDt(int finest_level, int, amrex::Vector<int> &n_cycle, const amrex::Vector<amrex::IntVect>&, amrex::Vector<amrex::Real> &dt_level, amrex::Real stop_time) override
Compute initial
dt
.- Parameters:
finest_level – Index of finest level
subcycle – unused
n_cycle – the number of subcycles at this level
ref_ratio – unused
dt_level – Real Vector, will contain initial timestep at each level.
stop_time – End time of simulation
-
void computeNewDt(int finest_level, int, amrex::Vector<int> &n_cycle, const amrex::Vector<amrex::IntVect>&, amrex::Vector<amrex::Real> &dt_min, amrex::Vector<amrex::Real> &dt_level, amrex::Real stop_time, int post_regrid_flag) override
Compute new
dt
.- Parameters:
finest_level – Index of finest level
subcycle – unused
n_cycle – the number of subcycles at this level
ref_ratio – unused
dt_min –
dt_level –
stop_time – End time of simulation
post_regrid_flag – Have we done regridding yet?
-
void allocOldData() override
Allocate data at old time.
-
void removeOldData() override
Remove data at old time.
-
void do_energy_diagnostics()
Print information about energy budget.
-
void post_timestep(int iteration) override
Do work after timestep().
- Parameters:
iteration – where we are in the current AMR subcycle
-
void postCoarseTimeStep(amrex::Real cumtime) override
Contains operations to be done only after a full coarse timestep.
- Parameters:
cumtime –
-
void post_regrid(int lbase, int new_finest) override
Do work after regrid().
- Parameters:
lbase –
new_finest –
-
void post_init(amrex::Real stop_time) override
Do work after init().
- Parameters:
stop_time – End time of simulation
-
void scf_relaxation()
Initialize a model using the self-consistent field method.
-
void do_hscf_solve()
Perform the Hachisu SCF method.
-
void post_grown_restart()
Do work after restart with
grown_factor
> 1
-
void errorEst(amrex::TagBoxArray &tags, int clearval, int tagval, amrex::Real time, int n_error_buf = 0, int ngrow = 0) override
Error estimation for regridding.
- Parameters:
tags – TagBoxArray of tags
clearval – integer
tagval – integer
time – current time
n_error_buf – integer
ngrow – Number of ghost cells
-
void apply_problem_tags(amrex::TagBoxArray &tags, amrex::Real time)
Apply the problem_tagging routine.
- Parameters:
tags – TagBoxArray of tags
time – current time
-
void apply_tagging_restrictions(amrex::TagBoxArray &tags, amrex::Real time)
Apply any tagging restrictions that must be satisfied by all problems.
- Parameters:
tags – TagBoxArray of tags
time – current time
-
std::unique_ptr<amrex::MultiFab> derive(const std::string &name, amrex::Real time, int ngrow) override
Returns a MultiFab containing the derived data for this level. If ngrow>0 the MultiFab is built on the appropriately grown BoxArray.
- Parameters:
name – Name of derived data
time – Current time
ngrow – Number of ghost cells
-
void derive(const std::string &name, amrex::Real time, amrex::MultiFab &mf, int dcomp) override
This version of derive() fills the dcomp’th component of mf with the derived quantity.
- Parameters:
name – Name of quantity to derive
time – current time
mf – MultiFab to store derived quantity in
dcomp – index of component of
mf
to fill with derived quantity
-
void final_radiation_call(amrex::MultiFab &S_new, int iteration, int ncycle)
- Parameters:
S_new –
iteration – where we are in the current AMR subcycle
ncycle – the number of subcycles at this level
-
void reset_internal_energy(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &State, int ng)
Reset the internal energy
- Parameters:
State – Current state (MultiFab)
ng – number of ghost cells
-
void check_div_B(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &state)
Check if divergence of B is zero at initialization
- Parameters:
Bx – magnetic field in x
By – magnetic field in y
Bz – magnetic field in z
state – the state to operate on
-
void computeTemp(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &state, amrex::Real time, int ng)
Compute the current temperature
- Parameters:
state – the state to operate on
time – current time
ng – number of ghost cells
-
void create_source_corrector()
Add any terms needed to correct the source terms. Currently this is a lagged predictor for CTU and the time-centering correction for simplified SDC.
-
void swap_state_time_levels(const amrex::Real dt)
This is a hack to make sure that we only ever have new data for certain state types that only ever need new time data; by doing a swap now, we’ll guarantee that allocOldData() does nothing. We do this because we never need the old data, so we don’t want to allocate memory for it.
- Parameters:
dt – timestep
-
amrex::Real volWgtSum(const MultiFab &mf, int comp, bool local = false, bool finemask = true)
Volume weighted sum of given quantity
- Parameters:
mf – MultiFab containing quantity
comp – integer index of quantity in MultiFab
local – boolean, is sum local (over each patch) or over entire MultiFab?
finemask – boolean, should we build a mask to exclude finer levels?
-
amrex::Real volWgtSum(const std::string &name, amrex::Real time, bool local = false, bool finemask = true)
Volume weighted sum of given quantity
- Parameters:
name – Name of quantity
time – current time
local – boolean, is sum local (over each patch) or over entire MultiFab?
finemask – boolean, should we build a mask to exclude finer levels?
-
amrex::Real locWgtSum(const MultiFab &mf, int comp, int idir, bool local = false)
Sum weighted by volume multiplied by distance from center in given direction
- Parameters:
mf – MultiFab containing quantity
comp – integer index of quantity in MultiFab
idir – Axis along which to compute distance from center
local – boolean, is sum local (over each patch) or over entire MultiFab?
-
amrex::Real locWgtSum(const std::string &name, amrex::Real time, int idir, bool local = false)
Sum weighted by volume multiplied by distance from center in given direction
- Parameters:
name – Name of quantity
time – current time
idir – Axis along which to compute distance from center
local – boolean, is sum local (over each patch) or over entire MultiFab?
-
amrex::Real volProductSum(const MultiFab &mf1, const MultiFab &mf2, int comp1, int comp2, bool local = false)
Volume weighted sum of the product of two quantities
- Parameters:
mf1 – MultiFab containing first quantity
mf2 – MultiFab containing second quantity
comp1 – integer index of quantity in first MultiFab
comp2 – integer index of quantity in second MultiFab
local – boolean, is sum local (over each patch) or over entire MultiFab?
-
amrex::Real volProductSum(const std::string &name1, const std::string &name2, amrex::Real time, bool local = false)
Volume weighted sum of the product of two quantities
- Parameters:
name1 – Name of first quantity
name2 – Name of second quantity
time – current time
local – boolean, is sum local (over each patch) or over entire MultiFab?
-
amrex::Real locSquaredSum(const std::string &name, amrex::Real time, int idir, bool local = false)
Location weighted sum of (quantity) squared
- Parameters:
name – Name of quantity
time – current time
idir – Direction along which to weight
local – boolean, is sum local (over each patch) or over entire MultiFab?
-
void gwstrain(Real time, Real &h_plus_1, Real &h_cross_1, Real &h_plus_2, Real &h_cross_2, Real &h_plus_3, Real &h_cross_3, bool local)
Calculate the gravitational wave signal
- Parameters:
time – current time
h_plus_1 – gravitational wave signal in x-direction
h_cross_1 – gravitational wave signal in x-direction
h_plus_2 – gravitational wave signal in y-direction
h_cross_2 – gravitational wave signal in y-direction
h_plus_3 – gravitational wave signal in z-direction
h_cross_3 – gravitational wave signal in z-direction
local – is sum local
-
void expand_state(amrex::MultiFab &S, amrex::Real time, int ng)
Fill a version of the state with
ng
ghost zones from the state data. After we do the fill, optionally clean the data using the iclean parameter.- Parameters:
S – MultiFab to be filled
time – current time
ng – number of ghost cells
-
void FluxRegCrseInit()
-
void FluxRegFineAdd()
-
void reflux(int crse_level, int fine_level, bool in_post_timestep)
Do refluxing
- Parameters:
crse_level –
fine_level –
in_post_timestep –
-
void enforce_min_density(amrex::MultiFab &state, int ng)
This routine sets the density in
state
to be larger than the density floor. Note that it will operate everywhere onstate
, including ghost zones.- Parameters:
state – state
ng – number of ghost cells
-
advance_status check_for_negative_density()
After a hydro advance, check to see whether any density went too small.
-
void clean_state(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &state, amrex::Real time, int ng)
Given
State_Type
state data, perform a number of cleaning steps to make sure the data is sensible.- Parameters:
state – State data
time – current time
ng – number of ghost cells
-
void avgDown()
Average new state from
level+1
down tolevel
-
void avgDown(int state_indx)
Average given component of new state from
level+1
down tolevel
- Parameters:
state_indx – Index of state to average down
-
void buildMetrics()
-
void initMFs()
Initialize the MultiFabs and flux registers that live as class members.
-
void sum_integrated_quantities()
Calculates volume weight sums of variables and prints to screen
-
void problem_diagnostics()
Problem-specific diagnostics (called by sum_integrated_quantities)
-
void write_info()
-
void define_new_center(const amrex::MultiFab &S, amrex::Real time)
Find the position of the “center” by interpolating from data at cell centers
- Parameters:
S – Current state
time – current time
-
void write_center()
-
void pointmass_update(amrex::Real time, amrex::Real dt)
- Parameters:
time – current time
dt – timestep
-
inline bool oversubscribing()
Returns true if we’re on the GPU and the total memory on this level oversubscribes GPU memory.
Public Members
Public Static Functions
-
static void writeBuildInfo()
Dump build info
-
static void variableSetUp()
Define data descriptors.
::variableSetUp is called in the constructor of Amr.cpp, so it should get called every time we start or restart a job
-
static void ErrorSetUp()
Define tagging functions.
-
static void variableCleanUp()
Cleanup data descriptors at end of run.
::variableSetUp is in
Castro_setup.cpp
variableCleanUp
is called once at the end of a simulation
-
static void read_particle_params()
Read particle-related inputs
-
static void check_for_nan(const amrex::MultiFab &state, int check_ghost = 0)
Check for NaNs in the given MultiFab
- Parameters:
state – MultiFab to check
check_ghost – do we check the ghost cells?
-
static void problem_post_simulation(amrex::Vector<std::unique_ptr<AmrLevel>> &amr_level)
Do work at the end of the simulation - before the last outputs
- Parameters:
amr_level –
-
static void reset_internal_energy(const amrex::Box &bx, amrex::Array4<amrex::Real> const Bx, amrex::Array4<amrex::Real> const By, amrex::Array4<amrex::Real> const Bz, amrex::Array4<amrex::Real> const u)
Reset the internal energy
- Parameters:
bx – Box to update
u – Current state (Fab)
-
static void add_magnetic_e(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &state)
Add magnetic contribution to energy density
- Parameters:
Bx – magnetic field in x
By – magnetic field in y
Bz – magnetic field in z
state – the state to operate on
-
static inline int using_point_mass()
Are we using point mass gravity?
-
static void extern_init()
-
static void read_params()
-
static void normalize_species(amrex::MultiFab &S_new, int ng)
Normalize species fractions so they sum to 1
- Parameters:
S_new – State
ng – number of ghost cells
-
static void enforce_consistent_e(amrex::MultiFab &Bx, amrex::MultiFab &By, amrex::MultiFab &Bz, amrex::MultiFab &S)
Enforces \f(\ \rho E = \rho e + \frac{1}{2} \rho (u^2 + v^2 + w^2) \f)
- Parameters:
S – State
-
static void enforce_speed_limit(amrex::MultiFab &state, int ng)
Ensure the magnitude of the velocity is not larger than
castro.speed_limit
in any zone.- Parameters:
state – state
ng – number of ghost cells
-
static inline void stopJob()
Public Static Attributes
-
static amrex::Real num_zones_advanced
A record of how many cells we have advanced throughout the simulation. This is saved as a real because we will be storing the number of zones advanced as a ratio with the number of zones on the coarse grid (to help prevent the number from getting too large), and that may be a non-integer number.
-
static Vector<std::unique_ptr<std::fstream>> data_logs
diagnostics
-
static Vector<std::unique_ptr<std::fstream>> problem_data_logs
-
static params_t params
runtime parameters
Protected Attributes
-
amrex::Vector<std::unique_ptr<amrex::iMultiFab>> ib_mask
Build a mask that ghost cells overlapping with interior cells in the same multifab are set to 0, whereas others are set to 1.
-
amrex::Vector<std::unique_ptr<amrex::StateData>> prev_state
State data to hold if we want to do a retry.
-
bool keep_prev_state = {}
Flag for indicating that we want to save prev_state until the reflux.
-
int iteration
Keep track of which AMR iteration we’re on.
-
int in_retry
-
int num_subcycles_taken
-
bool FillPatchedOldState_ok
-
int sub_iteration
subcycle capability
-
int sub_ncycle
-
int sdc_iteration
sdc
-
int current_sdc_node
Protected Static Attributes
-
static std::vector<std::string> burn_weight_names
-
static int SDC_NODES
-
static int do_cxx_prob_initialize
For initialization of the C++ values of the runtime parameters.
-
static bool signalStopJob
-
static int radius_grow
-
static std::vector<std::string> err_list_names
-
static std::vector<int> err_list_ng
-
static int num_err_list_default
-
static int NUM_GROW
-
static int NUM_GROW_SRC
-
static int lastDtPlotLimited
-
static amrex::Real previousCPUTimeUsed
for keeping track of the amount of CPU or GPU time used — this will persist after restarts
-
static int hydro_tile_size_has_been_tuned
-
static Long largest_box_from_hydro_tile_size_tuning
-
static int SDC_Source_Type
-
static int num_state_type
Friends
- friend class Radiation
-
Castro()