File Castro_hydro.H

Functions

void cons_to_prim(const amrex::Real time)

Calculate primitive variables from conserved variables (uses StateData)

Parameters:

time – current time

void cons_to_prim(amrex::MultiFab &u, amrex::MultiFab &q, amrex::MultiFab &qaux, const amrex::Real time)

Calculate primitive variables from given conserved variables

Parameters:
  • u – MultiFab of conserved variables

  • q – MultiFab to save primitive variables to

  • qaux – MultiFab of auxiliary variables

  • time – Current time

void cons_to_prim_fourth(const amrex::Real time)

convert the conservative state cell averages to primitive cell averages with 4th order accuracy

Parameters:

time – current time

void check_for_cfl_violation(const amrex::MultiFab &State, const amrex::Real dt)

Check to see if the CFL condition has been violated

Parameters:
  • State – MultiFab of conserved variables

  • dt – timestep

advance_status construct_ctu_hydro_source(amrex::Real time, amrex::Real dt)

this constructs the hydrodynamic source (essentially the flux divergence) using the CTU framework for unsplit hydrodynamics

Parameters:
  • time – current time

  • dt – timestep

void construct_mol_hydro_source(amrex::Real time, amrex::Real dt, amrex::MultiFab &A_update)

this constructs the hydrodynamic source (essentially the flux divergence) using method of lines integration. The output, is the divergence of the fluxes, A = -div{F(U)}

Parameters:
  • time – current time

  • dt – timestep

  • A_update – divergence of the fluxes

static void add_sdc_source_to_states(const Box &bx, const int idir, const Real dt, Array4<Real> const &qleft, Array4<Real> const &qright, Array4<Real const> const &sdc_src)

when using simplified-SDC, add the SDC source term predictor to the interface states

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction

  • dt – timestep

  • qleft – left state at the interface

  • qright – right state at the interface

  • sdc_source – primitive variable SDC source array

void src_to_prim(const amrex::Box &bx, const Real dt, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &old_src, amrex::Array4<amrex::Real const> const &src_corr, amrex::Array4<amrex::Real> const &srcQ)

convert the conserved form of the source terms, S(U), into source terms for the primitive variable equations, S(q).

Parameters:
  • bx – the box to operate over

  • q_arr – primitive variable state array

  • src – conserved source array

  • srcQ – primitive variable source array

static void ctoprim(const amrex::Box &bx, const amrex::Real time, amrex::Array4<amrex::Real const> const &uin, amrex::Array4<amrex::Real const> const &Bx, amrex::Array4<amrex::Real const> const &By, amrex::Array4<amrex::Real const> const &Bz, amrex::Array4<amrex::Real const> const &Erin, amrex::Array4<amrex::Real const> const &lam, amrex::Array4<amrex::Real> const &q_arr, amrex::Array4<amrex::Real> const &qaux_arr)

the actual work routine that does the conversion of conserved to primitive variables.

Parameters:
  • bx – the box to operate over

  • time – current time

  • uin – input conserved state

  • Bx – x-component of magentic field (if USE_MHD=TRUE)

  • By – y-component of magentic field (if USE_MHD=TRUE)

  • Bz – z-component of magentic field (if USE_MHD=TRUE)

  • Erin – radiation energy (if USE_RAD=TRUE)

  • lam – radiation flux limiter (if USE_RAD=TRUE)

  • q_arr – output primitive state

  • qaux_arr – output auxiliary quantities

void shock(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &U_src_arr, amrex::Array4<amrex::Real> const &shk)

A multidimensional shock detection algorithm

Parameters:
  • bx – the box to operate over

  • q_arr – the primitive variable state

  • U_scr_arr – the conservative state sources

  • shk – the shock flag (1 = shock, 0 = no shock)

void divu(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real> const &div)

Compute the node-centered velocity divergence (DU)_{i-1/2.j-1/2,k-1/2}

Parameters:
  • bx – the box to operate over

  • q_arr – input primitive variable state

  • div – output velocity divergence

void apply_av(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &div, amrex::Array4<amrex::Real const> const &uin, amrex::Array4<amrex::Real> const &flux)

Update the flux with the artificial viscosity. This is a 2nd order accurate implementation.

Parameters:
  • bx – the box to operate over

  • idir – the coordinate direction (0 = x, 1 = y, 2 = z)

  • div – the node-centered velocity divergence

  • uin – the conserved state

  • flux – the flux in direction idir

void apply_av_rad(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &div, amrex::Array4<amrex::Real const> const &Erin, amrex::Array4<amrex::Real> const &radflux)

Update the radiation flux with the artificial viscosity. This is a 2nd order accurate implementation.

Parameters:
  • bx – the box to operate over

  • idir – the coordinate direction (0 = x, 1 = y, 2 = z)

  • div – the node-centered velocity divergence

  • Erin – the radiation energy

  • radflux – the radiation energy flux in direction idir

void ctu_rad_consup(const amrex::Box &bx, amrex::Array4<amrex::Real> const &update, amrex::Array4<amrex::Real const> const &Erin, amrex::Array4<amrex::Real> const &Erout, amrex::Array4<amrex::Real const> const &radflux1, amrex::Array4<amrex::Real const> const &qx, amrex::Array4<amrex::Real const> const &area1, amrex::Array4<amrex::Real const> const &radflux2, amrex::Array4<amrex::Real const> const &qy, amrex::Array4<amrex::Real const> const &area2, amrex::Array4<amrex::Real const> const &radflux3, amrex::Array4<amrex::Real const> const &qz, amrex::Array4<amrex::Real const> const &area3, int &nstep_fsp, amrex::Array4<amrex::Real const> const &vol, const amrex::Real dt)
void ctu_ppm_states(const amrex::Box &bx, const amrex::Box &vbx, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real const> const &srcQ, amrex::Array4<amrex::Real> const &qxm, amrex::Array4<amrex::Real> const &qxp, amrex::Array4<amrex::Real> const &qym, amrex::Array4<amrex::Real> const &qyp, amrex::Array4<amrex::Real> const &qzm, amrex::Array4<amrex::Real> const &qzp, amrex::Array4<amrex::Real const> const &dloga, const amrex::Real dt)

Compute the normal left and right primitive variable interface states at each interface using the piecewise parabolic method with characteristic tracing. This implementation is for pure hydrodynamics

Parameters:
  • bx – the box to operate over

  • vbx – the valid region box (excludes ghost cells)

  • q_arr – the primitive variable state

  • qaux_arr – the auxiliary state

  • srcQ – the source terms for the primitive variable equations

  • qxm – left interface state in x, q_{i-1/2,j,k,L}

  • qxp – right interface state in x, q_{i-1/2,j,k,R}

  • qym – left interface state in y, q_{i,j-1/2,k,L}

  • qyp – right interface state in y, q_{i,j-1/2,k,R}

  • qzm – left interface state in z, q_{i,j,k-1/2,L}

  • qzp – right interface state in z, q_{i,j,k-1/2,R}

  • dloga – the geometric factor d(log area)

  • dt – timestep

void ctu_plm_states(const amrex::Box &bx, const amrex::Box &vbx, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real const> const &srcQ, amrex::Array4<amrex::Real> const &qxm, amrex::Array4<amrex::Real> const &qxp, amrex::Array4<amrex::Real> const &qym, amrex::Array4<amrex::Real> const &qyp, amrex::Array4<amrex::Real> const &qzm, amrex::Array4<amrex::Real> const &qzp, amrex::Array4<amrex::Real const> const &dloga, const amrex::Real dt)

Compute the normal left and right primitive variable interface states at each interface using piecewise linear reconstruction with characteristic tracing. This implementation is for pure hydrodynamics

Parameters:
  • bx – the box to operate over

  • vbx – the valid region box (excludes ghost cells)

  • q_arr – the primitive variable state

  • qaux_arr – the auxiliary state

  • srcQ – the source terms for the primitive variable equations

  • qxm – left interface state in x, q_{i-1/2,j,k,L}

  • qxp – right interface state in x, q_{i-1/2,j,k,R}

  • qym – left interface state in y, q_{i,j-1/2,k,L}

  • qyp – right interface state in y, q_{i,j-1/2,k,R}

  • qzm – left interface state in z, q_{i,j,k-1/2,L}

  • qzp – right interface state in z, q_{i,j,k-1/2,R}

  • dloga – the geometric factor d(log area)

  • dt – timestep

void ctu_ppm_rad_states(const amrex::Box &bx, const amrex::Box &vbx, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real const> const &srcQ, amrex::Array4<amrex::Real> const &qxm, amrex::Array4<amrex::Real> const &qxp, amrex::Array4<amrex::Real> const &qym, amrex::Array4<amrex::Real> const &qyp, amrex::Array4<amrex::Real> const &qzm, amrex::Array4<amrex::Real> const &qzp, amrex::Array4<amrex::Real const> const &dloga, const amrex::Real dt)

Compute the normal left and right primitive variable interface states at each interface using the piecewise parabolic method with characteristic tracing. This implementation is for the radiation hydrodynamics.

Parameters:
  • bx – the box to operate over

  • vbx – the valid region box (excludes ghost cells)

  • q_arr – the primitive variable state

  • qaux_arr – the auxiliary state

  • srcQ – the source terms for the primitive variable equations

  • qxm – left interface state in x, q_{i-1/2,j,k,L}

  • qxp – right interface state in x, q_{i-1/2,j,k,R}

  • qym – left interface state in y, q_{i,j-1/2,k,L}

  • qyp – right interface state in y, q_{i,j-1/2,k,R}

  • qzm – left interface state in z, q_{i,j,k-1/2,L}

  • qzp – right interface state in z, q_{i,j,k-1/2,R}

  • dloga – the geometric factor d(log area)

  • dt – timestep

void mol_plm_reconstruct(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &flatn_arr, amrex::Array4<amrex::Real const> const &src_q_arr, amrex::Array4<amrex::Real> const &dq, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp)

Compute the left and right primitive variable interface states at each interface using piecewise linear reconstruction for a method-of-lines type of integration (i.e., no prediction to n+1/2 or characteristic tracing). This implementation is for pure hydrodynamics

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction (0 = x, 1 = y, 2 = z)

  • q_arr – the primitive variable state

  • flatn_arr – flattening coefficient

  • dq – slope of the primitive variables

  • qm – left interface state, e.g., q_{i-1/2,j,k,L}

  • qp – right interface state, e.g., q_{i-1/2,j,k,R}

void mol_ppm_reconstruct(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &flatn_arr, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp)

Compute the left and right primitive variable interface states at each interface using piecewise parabolic reconstruction for a method-of-lines type of integration (i.e., no prediction to n+1/2 or characteristic tracing). This implementation is for pure hydrodynamics

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction (0 = x, 1 = y, 2 = z)

  • q_arr – the primitive variable state

  • flatn_arr – flattening coefficient

  • qm – left interface state, e.g., q_{i-1/2,j,k,L}

  • qp – right interface state, e.g., q_{i-1/2,j,k,R}

void mol_consup(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &shk, amrex::Array4<amrex::Real const> const &srcU, amrex::Array4<amrex::Real> const &update, const amrex::Real dt, amrex::Array4<amrex::Real const> const &flux0, amrex::Array4<amrex::Real const> const &flux1, amrex::Array4<amrex::Real const> const &flux2, amrex::Array4<amrex::Real const> const &area0, amrex::Array4<amrex::Real const> const &area1, amrex::Array4<amrex::Real const> const &area2, amrex::Array4<amrex::Real const> const &q0, amrex::Array4<amrex::Real const> const &q1, amrex::Array4<amrex::Real const> const &vol)

Compute the conservative update of the fluid state given the fluxes. This returns a source term, update, of the form -div{F} + S(U) which can then be added to the old state to update in time. This version is written for the method-of-lines / true SDC framework.

Parameters:
  • bx – the box to operate over

  • shk – the shock flag

  • srcU – source terms for the conservative equations

  • update – the final conserved state update, -div{F} + S(U)

  • flux0 – flux in the x-direction

  • flux1 – flux in the y-direction

  • flux2 – flux in the z-direction

  • area0 – area of x faces

  • area1 – area of y faces

  • area2 – area of z faces

  • q0 – Godunuv state on x faces

  • q1 – Godunuv state on y faces

  • vol – cell volume

void mol_diffusive_flux(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &U, amrex::Array4<amrex::Real const> const &cond, amrex::Array4<amrex::Real> const &flux)

Compute the diffusive flux term and add it to the energy fluxes. This implementation is second-order accurate.

Parameters:
  • bx – the box to operate over

  • idir – the coordinate direction (0 = x, 1 = y, 2 = z)

  • U – the conserved state

  • cond – thermal conductivity

  • flux – hydrodynamic flux in direction idir

void trans_single(const amrex::Box &bx, int idir_t, int idir_n, amrex::Array4<amrex::Real const> const &qm, amrex::Array4<amrex::Real> const &qmo, amrex::Array4<amrex::Real const> const &qp, amrex::Array4<amrex::Real> const &qpo, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real const> const &flux_t, amrex::Array4<amrex::Real const> const &rflux_t, amrex::Array4<amrex::Real const> const &q_t, amrex::Array4<amrex::Real const> const &area_t, amrex::Array4<amrex::Real const> const &vol, amrex::Real hdt, amrex::Real cdtdx)

Add the transverse flux difference in idir_t to the interface states in direction idir_n as part of the CTU hydrodynamics update. This is the main driver.

Parameters:
  • bx – the box to operate over

  • idir_t – direction for the transverse flux difference (0 = x, 1 = y, 2 = z)

  • idir_n – direction of the interface states normal (0 = x, 1 = y, 2 = z)

  • qm – input left interface state

  • qmo – updated left interface state

  • qp – input right interface state

  • qpo – updated right interface state

  • qaux – auxiliary state

  • flux_t – flux in the idir_t direction

  • rflux_t – radiation flux in the idir_t direction

  • q_t – Godunov state in the idir_t direction

  • area_t – face area in the idir_t direction

  • vol – cell volume

  • hdt – 1/2 * timestep

  • cdt – weight * timestep, where the weight comes from the CTU algorithm

void actual_trans_single(const amrex::Box &bx, int idir_t, int idir_n, int d, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real> const &qo_arr, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real const> const &flux_t, amrex::Array4<amrex::Real const> const &rflux_t, amrex::Array4<amrex::Real const> const &q_t, amrex::Array4<amrex::Real const> const &area_t, amrex::Array4<amrex::Real const> const &vol, amrex::Real hdt, amrex::Real cdtdx)

Add the transverse flux difference in idir_t to the interface states in direction idir_n as part of the CTU hydrodynamics update.

Parameters:
  • bx – the box to operate over

  • idir_t – direction for the transverse flux difference (0 = x, 1 = y, 2 = z)

  • idir_n – direction of the interface states normal (0 = x, 1 = y, 2 = z)

  • d – a flag set by the calling driver that determines which interface we are updating

  • q_arr – input interface state

  • qo_arr – updated interface state

  • qaux – auxiliary state

  • flux_t – flux in the idir_t direction

  • rflux_t – radiation flux in the idir_t direction

  • q_t – Godunov state in the idir_t direction

  • area_t – face area in the idir_t direction

  • vol – cell volume

  • hdt – 1/2 * timestep

  • cdt – weight * timestep, where the weight comes from the CTU algorithm

static void trans_final(const amrex::Box &bx, int idir_n, int idir_t1, int idir_t2, amrex::Array4<amrex::Real const> const &qm, amrex::Array4<amrex::Real> const &qmo, amrex::Array4<amrex::Real const> const &qp, amrex::Array4<amrex::Real> const &qpo, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real const> const &flux_t1, amrex::Array4<amrex::Real const> const &rflux_t1, amrex::Array4<amrex::Real const> const &flux_t2, amrex::Array4<amrex::Real const> const &rflux_t2, amrex::Array4<amrex::Real const> const &q_t1, amrex::Array4<amrex::Real const> const &q_t2, amrex::Real cdtdx_t1, amrex::Real cdtdx_t2)

The final transverse update where the flux difference in both perpendicular directions are added to the normal flux. This is the main driver.

Parameters:
  • bx – the box to operate over

  • idir_n – direction of the interface states normal (0 = x, 1 = y, 2 = z)

  • idir_t1 – direction for the first transverse flux difference (0 = x, 1 = y, 2 = z)

  • idir_t2 – direction for the second transverse flux difference (0 = x, 1 = y, 2 = z)

  • qm – input left interface state

  • qmo – updated left interface state

  • qp – input right interface state

  • qpo – updated right interface state

  • qaux – auxiliary state

  • flux_t1 – flux in the idir_t1 direction

  • rflux_t1 – radiation flux in the idir_t1 direction

  • flux_t2 – flux in the idir_t2 direction

  • rflux_t2 – radiation flux in the idir_t2 direction

  • q_t1 – Godunov state in the idir_t1 direction

  • q_t2 – Godunov state in the idir_t2 direction

  • hdt – 1/2 * timestep

  • cdtdx_t1 – weight * timestep in idir_t1 direction

  • cdtdx_t2 – weight * timestep in idir_t2 direction

static void actual_trans_final(const amrex::Box &bx, int idir_n, int idir_t1, int idir_t2, int d, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real> const &qo_arr, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real const> const &flux_t1, amrex::Array4<amrex::Real const> const &rflux_t1, amrex::Array4<amrex::Real const> const &flux_t2, amrex::Array4<amrex::Real const> const &rflux_t2, amrex::Array4<amrex::Real const> const &q_t1, amrex::Array4<amrex::Real const> const &q_t2, amrex::Real cdtdx_t1, amrex::Real cdtdx_t2)

The final transverse update where the flux difference in both perpendicular directions are added to the normal flux.

Parameters:
  • bx – the box to operate over

  • idir_n – direction of the interface states normal (0 = x, 1 = y, 2 = z)

  • idir_t1 – direction for the first transverse flux difference (0 = x, 1 = y, 2 = z)

  • idir_t2 – direction for the second transverse flux difference (0 = x, 1 = y, 2 = z)

  • d – a flag set by the driver to indicate which state we are operating on

  • q_arr – input interface state

  • qo_arr – updated interface state

  • qaux – auxiliary state

  • flux_t1 – flux in the idir_t1 direction

  • rflux_t1 – radiation flux in the idir_t1 direction

  • flux_t2 – flux in the idir_t2 direction

  • rflux_t2 – radiation flux in the idir_t2 direction

  • q_t1 – Godunov state in the idir_t1 direction

  • q_t2 – Godunov state in the idir_t2 direction

  • hdt – 1/2 * timestep

  • cdtdx_n – weight * timestep in idir_n direction

  • cdtdx_t1 – weight * timestep in idir_t1 direction

  • cdtdx_t2 – weight * timestep in idir_t2 direction

void trace_ppm(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real const> const &srcQ, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp, amrex::Array4<amrex::Real const> const &dloga, const amrex::Box &vbx, const amrex::Real dt)

Reconstruct the primitive state as parabola, integrate under them, and perform the characteristic tracing to get the interface states. This is for the CTU hydrodynamics scheme.

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction of the interface (0 = x, 1 = y, 2 = z)

  • q_arr – primitive variable state

  • qaux_arr – the auxiliary state

  • srcQ – primitive variable equation source terms

  • qm – left interface state (e.g., q_{i-1/2,j,k,L})

  • qp – right interface state (e.g., q_{i-1/2,j,k,R})

  • dloga – geometric factor dlog(area)

  • vbx – the valid region box (excluding ghost cells)

  • dt – timestep

void trace_plm(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp, amrex::Array4<amrex::Real const> const &dloga, amrex::Array4<amrex::Real const> const &SrcQ, const amrex::Box &vbx, const amrex::Real dt)

Reconstruct the primitive state as pieceeise linear, integrate under them, and perform the characteristic tracing to get the interface states. This is for the CTU hydrodynamics scheme.

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction of the interface (0 = x, 1 = y, 2 = z)

  • q_arr – primitive variable state

  • qaux_arr – the auxiliary state

  • qm – left interface state (e.g., q_{i-1/2,j,k,L})

  • qp – right interface state (e.g., q_{i-1/2,j,k,R})

  • dloga – geometric factor dlog(area)

  • srcQ – primitive variable equation source terms

  • vbx – the valid region box (excluding ghost cells)

  • dt – timestep

void trace_ppm_rad(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real const> const &U_arr, amrex::Array4<amrex::Real const> const &rho_inv_arr, amrex::Array4<amrex::Real const> const &q_arr, amrex::Array4<amrex::Real const> const &qaux_arr, amrex::Array4<amrex::Real const> const &srcQ, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp, amrex::Array4<amrex::Real const> const &dloga, const amrex::Box &vbx, const amrex::Real dt)

Reconstruct the primitive state as parabola, integrate under them, and perform the characteristic tracing to get the interface states. This is for the CTU radiation hydrodynamics scheme.

Parameters:
  • bx – the box to operate over

  • idir – coordinate direction of the interface (0 = x, 1 = y, 2 = z)

  • q_arr – primitive variable state

  • qaux_arr – the auxiliary state

  • srcQ – primitive variable equation source terms

  • qm – left interface state (e.g., q_{i-1/2,j,k,L})

  • qp – right interface state (e.g., q_{i-1/2,j,k,R})

  • dloga – geometric factor dlog(area)

  • vbx – the valid region box (excluding ghost cells)

  • dt – timestep

void consup_hydro(const amrex::Box &bx, amrex::Array4<amrex::Real> const &U_new, amrex::Array4<amrex::Real> const &flux0, amrex::Array4<amrex::Real const> const &qx, amrex::Array4<amrex::Real> const &flux1, amrex::Array4<amrex::Real const> const &qy, amrex::Array4<amrex::Real> const &flux2, amrex::Array4<amrex::Real const> const &qz, const amrex::Real dt)

Compute the conservative update of the hydrodynamics state and store it in update so it can be applied later to advance the state

Parameters:
  • bx – the box to operate over

  • shk – the shock flag

  • U_new – the new hydrodynamics conserved state

  • flux0 – flux in the x direction

  • qx – Godunov state in the x direction

  • flux1 – flux in the y direction

  • qy – Godunov state in the y direction

  • flux2 – flux in the z direction

  • qz – Godunov state in the z direction

  • dt – timestep

void cmpflx_plus_godunov(const amrex::Box &bx, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp, amrex::Array4<amrex::Real> const &flx, amrex::Array4<amrex::Real> const &rflx, amrex::Array4<amrex::Real> const &qgdnv, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real const> const &shk, const int idir, const bool store_full_state)

A general interface to the pure hydrodynamics and radiation Riemann solvers which takes the left and right interface states, qm and qp, and computes the flux through the interface and the godunov state (qint) on the interface

Parameters:
  • bx – the box to operate over

  • qm – left state on the interface

  • qp – right state on the interface

  • flx – flux through the interface

  • rflux – radiation fluxes through the interface

  • qgdnv – Godunov state on the interface (\either NQ or NGDNV)

  • qaux – auxiliary state

  • shk – shock flag

  • idir – coordinate direction of the solve (0 = x, 1 = y, 2 = z)

  • store_full_state – do we store all NQ or just the NGDNV subset in qgdnv

void compute_flux_from_q(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &qint, amrex::Array4<amrex::Real> const &F, const int idir)
static void store_godunov_state(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &qint, amrex::Array4<amrex::Real const> const &lambda, amrex::Array4<amrex::Real> const &qgdnv)
static void normalize_species_fluxes(const amrex::Box &bx, amrex::Array4<amrex::Real> const &flux)
static void limit_hydro_fluxes_on_small_dens(const amrex::Box &bx, int idir, amrex::Array4<amrex::Real const> const &u, amrex::Array4<amrex::Real const> const &vol, amrex::Array4<amrex::Real> const &flux, amrex::Array4<amrex::Real const> const &area, amrex::Real dt, bool scale_by_dAdt = true)
void enforce_reflect_states(const amrex::Box &bx, const int idir, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp)

For a reflecting boundary, we simply reflect the interface state just inside the domain to overwrite the state on the same interface but just outside the domain. This ensures that the flux through the reflecting boundary is zero.

void scale_flux(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &qint, amrex::Array4<amrex::Real> const &flux, amrex::Array4<amrex::Real const> const &area, const amrex::Real dt)
void scale_rad_flux(const amrex::Box &bx, amrex::Array4<amrex::Real> const &rflux, amrex::Array4<amrex::Real const> const &area, const amrex::Real dt)
static void reset_edge_state_thermo(const amrex::Box &bx, amrex::Array4<amrex::Real> const &qedge)
static void edge_state_temp_to_pres(const amrex::Box &bx, amrex::Array4<amrex::Real> const &qm, amrex::Array4<amrex::Real> const &qp)
void do_enforce_minimum_density(const amrex::Box &bx, amrex::Array4<amrex::Real> const &state_arr, const int verbose_warnings)
void construct_old_hybrid_source(amrex::MultiFab &source, amrex::MultiFab &state, amrex::Real time, amrex::Real dt)

Construct hybrid source terms at old time

Parameters:
  • source – MultiFab to save source to

  • state – Old state

  • time – current time

  • dt – timestep

void construct_new_hybrid_source(amrex::MultiFab &source, amrex::MultiFab &state_old, amrex::MultiFab &state_new, amrex::Real time, amrex::Real dt)

Construct hybrid source terms at new time

Parameters:
  • source – MultiFab to save source to

  • state_old – Old state

  • state_new – New state

  • time – current time

  • dt – timestep

void fill_hybrid_hydro_source(amrex::MultiFab &source, const amrex::MultiFab &state, const amrex::Real mult_factor)

Fill source with hybrid source terms

Parameters:
  • state – Current state

  • source – MultiFab to save source to

  • mult_factor

void hybrid_to_linear_momentum(amrex::MultiFab &state, int ng = 0)

Synchronize linear momentum with hybrid momentum

Parameters:
  • state – Current state

  • ng – Number of ghost cells

void linear_to_hybrid_momentum(amrex::MultiFab &state, int ng = 0)

Synchronize hybrid momentum with linear momentum

Parameters:
  • state – Current state

  • ng – Number of ghost cells

void fourth_interfaces(const amrex::Box &bx, const int idir, const int ncomp, amrex::Array4<amrex::Real const> const &a, amrex::Array4<amrex::Real> const &a_int)
void states(const amrex::Box &bx, const int idir, const int ncomp, amrex::Array4<amrex::Real const> const &a, amrex::Array4<amrex::Real const> const &a_int, amrex::Array4<amrex::Real const> const &flatn, amrex::Array4<amrex::Real> const &al, amrex::Array4<amrex::Real> const &ar)
void fourth_avisc(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &q, amrex::Array4<amrex::Real const> const &qaux, amrex::Array4<amrex::Real> const &avis, const int idir)
void fourth_add_diffusive_flux(const amrex::Box &bx, amrex::Array4<amrex::Real const> const &q, const int temp_comp, amrex::Array4<amrex::Real const> const &qint, amrex::Array4<amrex::Real> const &F, const int idir, const bool is_avg)
void make_cell_center(const Box &bx, Array4<Real const> const &U, Array4<Real> const &U_cc, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)
void make_cell_center_in_place(const Box &bx, Array4<Real> const &U, Array4<Real> const &tmp, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)
void compute_lap_term(const Box &bx, Array4<Real const> const &U, Array4<Real> const &lap, const int ncomp, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)
void make_fourth_average(const Box &bx, Array4<Real> const &q, Array4<Real const> const &q_bar, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)
void make_fourth_in_place(const Box &bx, Array4<Real> const &q, Array4<Real> const &tmp, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)
void make_fourth_in_place_n(const Box &bx, Array4<Real> const &q, const int ncomp, Array4<Real> const &tmp, GpuArray<int, 3> const &domlo, GpuArray<int, 3> const &domhi)