ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
WSM6 Class Reference

#include <ERF_WSM6.H>

Inheritance diagram for WSM6:
Collaboration diagram for WSM6:

Public Member Functions

 WSM6 ()
 
virtual ~WSM6 ()=default
 
void Define (SolverChoice &sc) override
 
void Init (const amrex::MultiFab &cons_in, const amrex::BoxArray &grids, const amrex::Geometry &geom, const amrex::Real &dt_advance, std::unique_ptr< amrex::MultiFab > &z_phys_nd, std::unique_ptr< amrex::MultiFab > &detJ_cc) override
 
void Set_dzmin (const amrex::Real dz_min) override
 
void Copy_State_to_Micro (const amrex::MultiFab &cons_in) override
 
void Copy_Micro_to_State (amrex::MultiFab &cons_in) override
 
void Update_Micro_Vars (amrex::MultiFab &cons_in) override
 
void Update_State_Vars (amrex::MultiFab &cons_in) override
 
void Advance (const amrex::Real &dt_advance, const SolverChoice &solverChoice) override
 
amrex::MultiFab * Qmoist_Ptr (const int &varIdx) override
 
int Qmoist_Size () override
 
int Qstate_Moist_Size () override
 
void Set_RealWidth (const int real_width) override
 
void Qmoist_Restart_Vars (const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
 
- Public Member Functions inherited from NullMoist
 NullMoist ()
 
virtual ~NullMoist ()=default
 
virtual int Qstate_NonMoist_Size ()
 
virtual void GetPlotVarNames (amrex::Vector< std::string > &a_vec) const
 
virtual void GetPlotVar (const std::string &, amrex::MultiFab &) const
 

Private Types

using FabPtr = std::shared_ptr< amrex::MultiFab >
 

Private Attributes

int m_qmoist_size = 3
 
int n_qstate_moist_size = 6
 
amrex::Vector< int > MicVarMap
 
amrex::Geometry m_geom
 
int m_real_width {0}
 
amrex::Real dt {0.0}
 
amrex::Real m_dzmin {0.0}
 
int nlev {0}
 
int zlo {0}
 
int zhi {0}
 
int m_axis {2}
 
bool m_do_cond {true}
 
amrex::MultiFab * m_z_phys_nd {nullptr}
 
amrex::MultiFab * m_detJ_cc {nullptr}
 
amrex::Array< FabPtr, MicVar_WSM6::NumVarsmic_fab_vars
 

Member Typedef Documentation

◆ FabPtr

using WSM6::FabPtr = std::shared_ptr<amrex::MultiFab>
private

Constructor & Destructor Documentation

◆ WSM6()

WSM6::WSM6 ( )
inline
41 {}

◆ ~WSM6()

virtual WSM6::~WSM6 ( )
virtualdefault

Member Function Documentation

◆ Advance()

void WSM6::Advance ( const amrex::Real dt_advance,
const SolverChoice solverChoice 
)
overridevirtual

Reimplemented from NullMoist.

10 {
11  dt = dt_advance;
12 
13 #ifdef ERF_USE_WSM6_FORT
14  static bool wsm6_inited = false;
15 
16  // Minimal phase-1 initialization for single-moment WSM6.
17  if (!wsm6_inited) {
18  constexpr double den0 = 1.28; // Standard dry-air density (kg/m^3)
19  constexpr double denr = static_cast<double>(rhoh2o);
20  constexpr double dens = static_cast<double>(rhos);
21  constexpr double cl = static_cast<double>(Cp_l);
22  constexpr double cpv = static_cast<double>(Cp_v);
23  constexpr int hail_opt = 0; // Graupel mode
24  mp_wsm6_init_c(den0, denr, dens, cl, cpv, hail_opt);
25  wsm6_inited = true;
26  }
27 
28  constexpr double g = static_cast<double>(CONST_GRAV);
29  constexpr double cpd = static_cast<double>(Cp_d);
30  constexpr double cpv = static_cast<double>(Cp_v);
31  constexpr double rd = static_cast<double>(R_d);
32  constexpr double rv = static_cast<double>(R_v);
33  constexpr double t0c = 273.15;
34  constexpr double ep1 = static_cast<double>(R_v / R_d - one);
35  constexpr double ep2 = static_cast<double>(R_d / R_v);
36  constexpr double qmin = 1.0e-12;
37  constexpr double xls = static_cast<double>(lsub);
38  constexpr double xlv0 = static_cast<double>(lat_vap);
39  constexpr double xlf0 = static_cast<double>(lat_ice);
40  constexpr double den0 = 1.28;
41  constexpr double denr = static_cast<double>(rhoh2o);
42  constexpr double cliq = static_cast<double>(Cp_l);
43  constexpr double cice = 2106.0;
44  constexpr double psat = 610.78;
45 
46  for (MFIter mfi(*mic_fab_vars[MicVar_WSM6::qv], TileNoZ()); mfi.isValid(); ++mfi) {
47  const Box box = mfi.tilebox();
48  const Box fab_box = mfi.fabbox();
49 
50  auto const& t_arr = mic_fab_vars[MicVar_WSM6::tabs]->array(mfi);
51  auto const& qv_arr = mic_fab_vars[MicVar_WSM6::qv]->array(mfi);
52  auto const& qc_arr = mic_fab_vars[MicVar_WSM6::qc]->array(mfi);
53  auto const& qi_arr = mic_fab_vars[MicVar_WSM6::qi]->array(mfi);
54  auto const& qr_arr = mic_fab_vars[MicVar_WSM6::qr]->array(mfi);
55  auto const& qs_arr = mic_fab_vars[MicVar_WSM6::qs]->array(mfi);
56  auto const& qg_arr = mic_fab_vars[MicVar_WSM6::qg]->array(mfi);
57  auto const& den_arr = mic_fab_vars[MicVar_WSM6::rho]->array(mfi);
58  auto const& p_arr = mic_fab_vars[MicVar_WSM6::pres]->array(mfi);
59  auto rain_arr = mic_fab_vars[MicVar_WSM6::rain_accum]->array(mfi);
60  auto snow_arr = mic_fab_vars[MicVar_WSM6::snow_accum]->array(mfi);
61  auto graup_arr = mic_fab_vars[MicVar_WSM6::graup_accum]->array(mfi);
62 
63  const int ilo = box.smallEnd(0);
64  const int ihi = box.bigEnd(0);
65  const int jlo = box.smallEnd(1);
66  const int jhi = box.bigEnd(1);
67  const int klo = box.smallEnd(2);
68  const int khi = box.bigEnd(2);
69 
70  const int imlo = fab_box.smallEnd(0);
71  const int imhi = fab_box.bigEnd(0);
72  const int jmlo = fab_box.smallEnd(1);
73  const int jmhi = fab_box.bigEnd(1);
74  const int kmlo = fab_box.smallEnd(2);
75  const int kmhi = fab_box.bigEnd(2);
76 
77  const Real dz_val = m_geom.CellSize(m_axis);
78  FArrayBox delz_fab(fab_box, 1);
79  auto const& delz_arr = delz_fab.array();
80  delz_fab.setVal(dz_val);
81 
82  const Array4<const Real> z_arr = (m_z_phys_nd) ? m_z_phys_nd->const_array(mfi) : Array4<const Real> {};
83  ParallelFor(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
84  delz_arr(i,j,k) = (z_arr) ? Real(0.25) * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
85  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
86  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
87  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz_val;
88  });
89 
90  Box box2d(fab_box);
91  box2d.makeSlab(2, 0);
92  FArrayBox rainncv_fab(box2d, 1);
93  FArrayBox sr_fab(box2d, 1);
94  FArrayBox snowncv_fab(box2d, 1);
95  FArrayBox graupelncv_fab(box2d, 1);
96  FArrayBox rainacc_fab(box2d, 1);
97  FArrayBox snowacc_fab(box2d, 1);
98  FArrayBox graupacc_fab(box2d, 1);
99 
100  auto const& rainncv_arr = rainncv_fab.array();
101  auto const& sr_arr = sr_fab.array();
102  auto const& snowncv_arr = snowncv_fab.array();
103  auto const& graupelncv_arr = graupelncv_fab.array();
104  auto const& rainacc_arr = rainacc_fab.array();
105  auto const& snowacc_arr = snowacc_fab.array();
106  auto const& graupacc_arr = graupacc_fab.array();
107  ParallelFor(box2d, [=] AMREX_GPU_DEVICE (int i, int j, int) {
108  rainacc_arr(i,j,0) = rain_arr(i,j,klo);
109  snowacc_arr(i,j,0) = snow_arr(i,j,klo);
110  graupacc_arr(i,j,0) = graup_arr(i,j,klo);
111  rainncv_arr(i,j,0) = Real(0.0);
112  sr_arr(i,j,0) = Real(0.0);
113  snowncv_arr(i,j,0) = Real(0.0);
114  graupelncv_arr(i,j,0) = Real(0.0);
115  });
116 
118  t_arr.dataPtr(),
119  qv_arr.dataPtr(), qc_arr.dataPtr(), qi_arr.dataPtr(),
120  qr_arr.dataPtr(), qs_arr.dataPtr(), qg_arr.dataPtr(),
121  den_arr.dataPtr(), p_arr.dataPtr(), delz_arr.dataPtr(),
122  static_cast<double>(dt), g, cpd, cpv, rd, rv, t0c, ep1, ep2, qmin,
123  xls, xlv0, xlf0, den0, denr, cliq, cice, psat,
124  rainacc_arr.dataPtr(), rainncv_arr.dataPtr(), sr_arr.dataPtr(),
125  snowacc_arr.dataPtr(), snowncv_arr.dataPtr(),
126  graupacc_arr.dataPtr(), graupelncv_arr.dataPtr(),
127  imlo, imhi, jmlo, jmhi, kmlo, kmhi,
128  ilo, ihi, jlo, jhi, klo, khi);
129 
130  ParallelFor(box2d, [=] AMREX_GPU_DEVICE (int i, int j, int) {
131  rain_arr(i,j,klo) = rainacc_arr(i,j,0);
132  snow_arr(i,j,klo) = snowacc_arr(i,j,0);
133  graup_arr(i,j,klo) = graupacc_arr(i,j,0);
134  });
135  }
136 #else
137  amrex::Abort("WSM6 Fortran bridge requested but ERF was not built with ERF_USE_WSM6_FORT");
138 #endif
139 }
constexpr amrex::Real R_v
Definition: ERF_Constants.H:21
constexpr amrex::Real lat_vap
Definition: ERF_Constants.H:95
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:22
constexpr amrex::Real lsub
Definition: ERF_Constants.H:78
constexpr amrex::Real rhoh2o
Definition: ERF_Constants.H:104
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real lat_ice
Definition: ERF_Constants.H:96
constexpr amrex::Real rhos
Definition: ERF_Constants.H:39
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:31
constexpr amrex::Real Cp_l
Definition: ERF_Constants.H:24
constexpr amrex::Real Cp_v
Definition: ERF_Constants.H:23
constexpr amrex::Real R_d
Definition: ERF_Constants.H:20
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
auto qv_arr
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:210
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
void mp_wsm6_run_c(double *t, double *qv, double *qc, double *qi, double *qr, double *qs, double *qg, double *den, double *p, double *delz, double delt, double g, double cpd, double cpv, double rd, double rv, double t0c, double ep1, double ep2, double qmin, double xls, double xlv0, double xlf0, double den0, double denr, double cliq, double cice, double psat, double *rain, double *rainncv, double *sr, double *snow, double *snowncv, double *graupel, double *graupelncv, int ims, int ime, int jms, int jme, int kms, int kme, int its, int ite, int jts, int jte, int kts, int kte)
void mp_wsm6_init_c(double den0, double denr, double dens, double cl, double cpv, int hail_opt)
amrex::MultiFab * m_z_phys_nd
Definition: ERF_WSM6.H:108
int m_axis
Definition: ERF_WSM6.H:105
amrex::Geometry m_geom
Definition: ERF_WSM6.H:100
amrex::Real dt
Definition: ERF_WSM6.H:102
amrex::Array< FabPtr, MicVar_WSM6::NumVars > mic_fab_vars
Definition: ERF_WSM6.H:111
@ qr
Definition: ERF_WSM6.H:27
@ qi
Definition: ERF_WSM6.H:26
@ qs
Definition: ERF_WSM6.H:28
@ qv
Definition: ERF_WSM6.H:24
@ qc
Definition: ERF_WSM6.H:25
@ rain_accum
Definition: ERF_WSM6.H:30
@ snow_accum
Definition: ERF_WSM6.H:31
@ rho
Definition: ERF_WSM6.H:20
@ graup_accum
Definition: ERF_WSM6.H:32
@ qg
Definition: ERF_WSM6.H:29
@ pres
Definition: ERF_WSM6.H:23
@ tabs
Definition: ERF_WSM6.H:22
real(c_double), parameter cice
Definition: ERF_module_model_constants.F90:30
real(c_double), parameter g
Definition: ERF_module_model_constants.F90:19
real(c_double), parameter cliq
Definition: ERF_module_model_constants.F90:29
real(c_double), parameter xlv0
Definition: ERF_module_model_constants.F90:51
real(c_double), parameter cpv
Definition: ERF_module_model_constants.F90:26
real(c_double), parameter xls
Definition: ERF_module_model_constants.F90:56
real(c_double), parameter psat
Definition: ERF_module_model_constants.F90:31
real(kind=kind_phys), parameter, private dens
Definition: ERF_module_mp_wsm6.F90:39
Here is the call graph for this function:

◆ Copy_Micro_to_State()

void WSM6::Copy_Micro_to_State ( amrex::MultiFab &  cons_in)
overridevirtual

Reimplemented from NullMoist.

7 {
8  for (MFIter mfi(cons, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
9  const auto& box3d = mfi.tilebox();
10  auto states = cons.array(mfi);
11 
12  auto rho = mic_fab_vars[MicVar_WSM6::rho]->array(mfi);
13  auto theta = mic_fab_vars[MicVar_WSM6::theta]->array(mfi);
14 
15  auto qv = mic_fab_vars[MicVar_WSM6::qv]->array(mfi);
16  auto qc = mic_fab_vars[MicVar_WSM6::qc]->array(mfi);
17  auto qi = mic_fab_vars[MicVar_WSM6::qi]->array(mfi);
18  auto qr = mic_fab_vars[MicVar_WSM6::qr]->array(mfi);
19  auto qs = mic_fab_vars[MicVar_WSM6::qs]->array(mfi);
20  auto qg = mic_fab_vars[MicVar_WSM6::qg]->array(mfi);
21 
22  ParallelFor(box3d, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
23  states(i,j,k,RhoTheta_comp) = rho(i,j,k) * theta(i,j,k);
24  states(i,j,k,RhoQ1_comp) = rho(i,j,k) * amrex::max(Real(0), qv(i,j,k));
25  states(i,j,k,RhoQ2_comp) = rho(i,j,k) * amrex::max(Real(0), qc(i,j,k));
26  states(i,j,k,RhoQ3_comp) = rho(i,j,k) * amrex::max(Real(0), qi(i,j,k));
27  states(i,j,k,RhoQ4_comp) = rho(i,j,k) * amrex::max(Real(0), qr(i,j,k));
28  states(i,j,k,RhoQ5_comp) = rho(i,j,k) * amrex::max(Real(0), qs(i,j,k));
29  states(i,j,k,RhoQ6_comp) = rho(i,j,k) * amrex::max(Real(0), qg(i,j,k));
30  });
31  }
32 
33  cons.FillBoundary(m_geom.periodicity());
34 }
#define RhoQ4_comp
Definition: ERF_IndexDefines.H:45
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ3_comp
Definition: ERF_IndexDefines.H:44
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoQ6_comp
Definition: ERF_IndexDefines.H:47
#define RhoQ5_comp
Definition: ERF_IndexDefines.H:46
rho
Definition: ERF_InitCustomPert_Bubble.H:106
@ theta
Definition: ERF_MM5.H:20
@ qv
Definition: ERF_Kessler.H:28
@ qc
Definition: ERF_SatAdj.H:36
@ theta
Definition: ERF_WSM6.H:21
@ cons
Definition: ERF_IndexDefines.H:158

Referenced by Update_State_Vars().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Copy_State_to_Micro()

void WSM6::Copy_State_to_Micro ( const amrex::MultiFab &  cons_in)
overridevirtual

Reimplemented from NullMoist.

38 {
39  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
40  // Match Morrison behavior: refresh microphysics ghost zones from state.
41  // WSM6 Fortran reads the full (ims:ime, jms:jme, kms:kme) slab.
42  const auto& box3d = mfi.growntilebox();
43  auto states = cons_in.array(mfi);
44 
45  auto rho = mic_fab_vars[MicVar_WSM6::rho]->array(mfi);
46  auto theta = mic_fab_vars[MicVar_WSM6::theta]->array(mfi);
47  auto tabs = mic_fab_vars[MicVar_WSM6::tabs]->array(mfi);
48  auto pres = mic_fab_vars[MicVar_WSM6::pres]->array(mfi);
49 
50  auto qv = mic_fab_vars[MicVar_WSM6::qv]->array(mfi);
51  auto qc = mic_fab_vars[MicVar_WSM6::qc]->array(mfi);
52  auto qi = mic_fab_vars[MicVar_WSM6::qi]->array(mfi);
53  auto qr = mic_fab_vars[MicVar_WSM6::qr]->array(mfi);
54  auto qs = mic_fab_vars[MicVar_WSM6::qs]->array(mfi);
55  auto qg = mic_fab_vars[MicVar_WSM6::qg]->array(mfi);
56 
57  ParallelFor(box3d, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
58  rho(i,j,k) = states(i,j,k,Rho_comp);
59  theta(i,j,k) = states(i,j,k,RhoTheta_comp) / states(i,j,k,Rho_comp);
60 
61  qv(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ1_comp) / states(i,j,k,Rho_comp));
62  qc(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ2_comp) / states(i,j,k,Rho_comp));
63  qi(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ3_comp) / states(i,j,k,Rho_comp));
64  qr(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ4_comp) / states(i,j,k,Rho_comp));
65  qs(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ5_comp) / states(i,j,k,Rho_comp));
66  qg(i,j,k) = amrex::max(Real(0.0), states(i,j,k,RhoQ6_comp) / states(i,j,k,Rho_comp));
67 
68  tabs(i,j,k) = getTgivenRandRTh(states(i,j,k,Rho_comp),
69  states(i,j,k,RhoTheta_comp),
70  qv(i,j,k));
71  pres(i,j,k) = getPgivenRTh(states(i,j,k,RhoTheta_comp), qv(i,j,k));
72  });
73  }
74 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
#define Rho_comp
Definition: ERF_IndexDefines.H:36
@ tabs
Definition: ERF_Kessler.H:24
@ pres
Definition: ERF_Kessler.H:25

Referenced by Update_Micro_Vars().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Define()

void WSM6::Define ( SolverChoice sc)
inlineoverridevirtual

Reimplemented from NullMoist.

45  {
46  m_axis = sc.ave_plane;
47  m_do_cond = (!sc.use_shoc);
48  }
bool m_do_cond
Definition: ERF_WSM6.H:106
bool use_shoc
Definition: ERF_DataStruct.H:1187
int ave_plane
Definition: ERF_DataStruct.H:1220

◆ Init()

void WSM6::Init ( const amrex::MultiFab &  cons_in,
const amrex::BoxArray &  grids,
const amrex::Geometry &  geom,
const amrex::Real dt_advance,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
std::unique_ptr< amrex::MultiFab > &  detJ_cc 
)
overridevirtual

Reimplemented from NullMoist.

15 {
16  dt = dt_advance;
17  m_geom = geom;
18 
19  m_z_phys_nd = z_phys_nd.get();
20  m_detJ_cc = detJ_cc.get();
21 
22  MicVarMap.resize(m_qmoist_size);
24 
25  for (int ivar = 0; ivar < MicVar_WSM6::NumVars; ++ivar) {
26  mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
27  1, cons_in.nGrowVect());
28  mic_fab_vars[ivar]->setVal(0.0);
29  }
30 
31  nlev = m_geom.Domain().length(2);
32  zlo = m_geom.Domain().smallEnd(2);
33  zhi = m_geom.Domain().bigEnd(2);
34 }
amrex::MultiFab * m_detJ_cc
Definition: ERF_WSM6.H:109
int zlo
Definition: ERF_WSM6.H:104
int nlev
Definition: ERF_WSM6.H:104
int zhi
Definition: ERF_WSM6.H:104
int m_qmoist_size
Definition: ERF_WSM6.H:95
amrex::Vector< int > MicVarMap
Definition: ERF_WSM6.H:98
@ NumVars
Definition: ERF_WSM6.H:33

◆ Qmoist_Ptr()

amrex::MultiFab* WSM6::Qmoist_Ptr ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullMoist.

76  {
78  return mic_fab_vars[MicVarMap[varIdx]].get();
79  }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
Here is the call graph for this function:

◆ Qmoist_Restart_Vars()

void WSM6::Qmoist_Restart_Vars ( const SolverChoice ,
std::vector< int > &  a_idx,
std::vector< std::string > &  a_names 
) const
inlineoverridevirtual

Reimplemented from NullMoist.

89  {
90  a_idx = {0, 1, 2};
91  a_names = {"RainAccum", "SnowAccum", "GraupAccum"};
92  }

◆ Qmoist_Size()

int WSM6::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

81 { return m_qmoist_size; }

◆ Qstate_Moist_Size()

int WSM6::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

82 { return n_qstate_moist_size; }
int n_qstate_moist_size
Definition: ERF_WSM6.H:96

◆ Set_dzmin()

void WSM6::Set_dzmin ( const amrex::Real  dz_min)
inlineoverridevirtual

Reimplemented from NullMoist.

57 { m_dzmin = dz_min; }
amrex::Real m_dzmin
Definition: ERF_WSM6.H:103

◆ Set_RealWidth()

void WSM6::Set_RealWidth ( const int  real_width)
inlineoverridevirtual

Reimplemented from NullMoist.

84 { m_real_width = real_width; }
int m_real_width
Definition: ERF_WSM6.H:101

◆ Update_Micro_Vars()

void WSM6::Update_Micro_Vars ( amrex::MultiFab &  cons_in)
inlineoverridevirtual

Reimplemented from NullMoist.

63  {
64  Copy_State_to_Micro(cons_in);
65  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitWSM6.cpp:37
Here is the call graph for this function:

◆ Update_State_Vars()

void WSM6::Update_State_Vars ( amrex::MultiFab &  cons_in)
inlineoverridevirtual

Reimplemented from NullMoist.

68  {
69  Copy_Micro_to_State(cons_in);
70  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateWSM6.cpp:6
Here is the call graph for this function:

Member Data Documentation

◆ dt

amrex::Real WSM6::dt {0.0}
private

◆ m_axis

int WSM6::m_axis {2}
private

Referenced by Define().

◆ m_detJ_cc

amrex::MultiFab* WSM6::m_detJ_cc {nullptr}
private

◆ m_do_cond

bool WSM6::m_do_cond {true}
private

Referenced by Define().

◆ m_dzmin

amrex::Real WSM6::m_dzmin {0.0}
private

Referenced by Set_dzmin().

◆ m_geom

amrex::Geometry WSM6::m_geom
private

◆ m_qmoist_size

int WSM6::m_qmoist_size = 3
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_real_width

int WSM6::m_real_width {0}
private

Referenced by Set_RealWidth().

◆ m_z_phys_nd

amrex::MultiFab* WSM6::m_z_phys_nd {nullptr}
private

◆ mic_fab_vars

amrex::Array<FabPtr, MicVar_WSM6::NumVars> WSM6::mic_fab_vars
private

Referenced by Qmoist_Ptr().

◆ MicVarMap

amrex::Vector<int> WSM6::MicVarMap
private

Referenced by Qmoist_Ptr().

◆ n_qstate_moist_size

int WSM6::n_qstate_moist_size = 6
private

Referenced by Qstate_Moist_Size().

◆ nlev

int WSM6::nlev {0}
private

◆ zhi

int WSM6::zhi {0}
private

◆ zlo

int WSM6::zlo {0}
private

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