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

#include <ERF_SatAdj.H>

Inheritance diagram for SatAdj:
Collaboration diagram for SatAdj:

Public Member Functions

 SatAdj ()
 
virtual ~SatAdj ()=default
 
void AdvanceSatAdj (const SolverChoice &)
 
void Define (SolverChoice &sc) override
 
void Init (const amrex::MultiFab &cons_in, const amrex::BoxArray &, const amrex::Geometry &geom, const amrex::Real &dt_advance, std::unique_ptr< amrex::MultiFab > &, std::unique_ptr< amrex::MultiFab > &) 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, const amrex::MultiFab &) 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 Qmoist_Restart_Vars (const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
 
void Set_RealWidth (const int real_width) 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
 
virtual void GetPlotVar (const std::string &a_name, amrex::MultiFab &a_mf, const int) const
 
virtual void SetCurrentLevel (const int &)
 
virtual void InitLevel (const int, const amrex::MultiFab &)
 
virtual int getDiagnosticsInterval () const
 
virtual void Set_dzmin (const amrex::Real)
 

Static Public Member Functions

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonSolveSatTemperature (const amrex::Real fac_cond, const amrex::Real tabs_old, const amrex::Real pres_mbar, const amrex::Real qv_old, amrex::Real &qsat)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void AdjustSatAdjCell (const amrex::Real fac_cond, const amrex::Real rdOcp, amrex::Real &tabs, const amrex::Real pres_mbar, amrex::Real &theta, amrex::Real &qv, amrex::Real &qc)
 

Private Types

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

Private Attributes

int m_qmoist_size = 0
 
int n_qstate_moist_size = 2
 
amrex::Geometry m_geom
 
int m_real_width {0}
 
amrex::Real dt
 
amrex::Real m_fac_cond
 
amrex::Real m_rdOcp
 
bool m_do_cond
 
amrex::Array< FabPtr, MicVar_SatAdj::NumVarsmic_fab_vars
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ SatAdj()

SatAdj::SatAdj ( )
inline
47 {}

◆ ~SatAdj()

virtual SatAdj::~SatAdj ( )
virtualdefault

Member Function Documentation

◆ AdjustSatAdjCell()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void SatAdj::AdjustSatAdjCell ( const amrex::Real  fac_cond,
const amrex::Real  rdOcp,
amrex::Real tabs,
const amrex::Real  pres_mbar,
amrex::Real theta,
amrex::Real qv,
amrex::Real qc 
)
inlinestatic
190  {
191  amrex::Real qsat;
192  erf_qsatw(tabs, pres_mbar, qsat);
193 
194  if ((qv + qc) > qsat) {
195 #if defined(AMREX_DEBUG)
196  const amrex::Real qvprev = qv;
197  const amrex::Real qcprev = qc;
198 #endif
199 
200  if (qc < 0) {
201  qv += qc;
202  qc = zero;
203  }
204 
205  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
206 
207  const amrex::Real delta_qv = qv - qsat;
208  qv = qsat;
209  qc += delta_qv;
210 
211 #if defined(AMREX_DEBUG)
212  amrex::Real qsatnew;
213  erf_qsatw(tabs, pres_mbar, qsatnew);
214  AMREX_ASSERT(std::abs(qv-qsatnew) < 1e-12);
215  AMREX_ASSERT(std::abs(qv+qc-qvprev-qcprev) < 1e-14);
216 #endif
217 
218  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
219  } else {
220  const amrex::Real delta_qc = qc;
221 
222  qv += qc;
223  qc = zero;
224 
225  tabs -= fac_cond * delta_qc;
226  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
227 
228  erf_qsatw(tabs, pres_mbar, qsat);
229  if (qv > qsat) {
230 #if defined(AMREX_DEBUG)
231  const amrex::Real qvprev = qv;
232  const amrex::Real qcprev = qc;
233  const amrex::Real Tprev = tabs;
234 #endif
235 
236  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
237 
238  const amrex::Real delta_qv = qv - qsat;
239  qv = qsat;
240  qc += delta_qv;
241 
242 #if defined(AMREX_DEBUG)
243  amrex::Real qsatnew;
244  erf_qsatw(tabs, pres_mbar, qsatnew);
245  AMREX_ASSERT(qv < qvprev);
246  AMREX_ASSERT(qc > qcprev);
247  AMREX_ASSERT(tabs > Tprev);
248  AMREX_ASSERT(std::abs(qv-qsatnew) < 1e-14);
249  AMREX_ASSERT(std::abs(qv+qc-qvprev-qcprev) < 1e-14);
250 #endif
251 
252  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
253  }
254  }
255  }
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
const Real rdOcp
Definition: ERF_InitCustomPert_Bomex.H:16
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:228
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonSolveSatTemperature(const amrex::Real fac_cond, const amrex::Real tabs_old, const amrex::Real pres_mbar, const amrex::Real qv_old, amrex::Real &qsat)
Definition: ERF_SatAdj.H:134
@ theta
Definition: ERF_SatAdj.H:31
@ qv
Definition: ERF_SatAdj.H:35
@ tabs
Definition: ERF_SatAdj.H:32
@ qc
Definition: ERF_SatAdj.H:36
Here is the call graph for this function:

◆ Advance()

void SatAdj::Advance ( const amrex::Real dt_advance,
const SolverChoice solverChoice 
)
inlineoverridevirtual

Reimplemented from NullMoist.

100  {
101  dt = dt_advance;
102 
103  this->AdvanceSatAdj(solverChoice);
104  }
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:8
amrex::Real dt
Definition: ERF_SatAdj.H:271
Here is the call graph for this function:

◆ AdvanceSatAdj()

void SatAdj::AdvanceSatAdj ( const SolverChoice )

Compute Precipitation-related Microphysics quantities.

9 {
10  if (!m_do_cond) { return; }
11 
13 
14  // Expose for GPU
15  Real d_fac_cond = m_fac_cond;
16  Real rdOcp = m_rdOcp;
17 
18  auto domain = m_geom.Domain();
19  int i_lo = domain.smallEnd(0);
20  int i_hi = domain.bigEnd(0);
21  int j_lo = domain.smallEnd(1);
22  int j_hi = domain.bigEnd(1);
23 
24  // get the temperature, density, theta, qt and qc from input
25  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
26 
27  auto tbx = mfi.tilebox();
28  if (tbx.smallEnd(0) == i_lo) { tbx.growLo(0,-m_real_width); }
29  if (tbx.bigEnd(0) == i_hi) { tbx.growHi(0,-m_real_width); }
30  if (tbx.smallEnd(1) == j_lo) { tbx.growLo(1,-m_real_width); }
31  if (tbx.bigEnd(1) == j_hi) { tbx.growHi(1,-m_real_width); }
32 
33  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
34  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
35  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
36  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
37  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
38 
39  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
40  {
41  Real T = tabs_array(i,j,k);
42  Real p = pres_array(i,j,k);
43  Real th = theta_array(i,j,k);
44  Real qv = qv_array(i,j,k);
45  Real qc = qc_array(i,j,k);
46 
47  AdjustSatAdjCell(d_fac_cond, rdOcp, T, p, th, qv, qc);
48 
49  tabs_array(i,j,k) = T;
50  theta_array(i,j,k) = th;
51  qv_array(i,j,k) = qv;
52  qc_array(i,j,k) = qc;
53  });
54  }
55 }
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+myhalf) *dx[0];const Real y=(j+myhalf) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - one)/(two *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(one+deltaT, inv_gm1);const Real T=(one+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=fourth *(one+std::cos(PI *std::min(r2d_xy, R)/R));})
Real * p
Definition: ERF_InitCustomPert_SquallLine.H:61
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE void AdjustSatAdjCell(const amrex::Real fac_cond, const amrex::Real rdOcp, amrex::Real &tabs, const amrex::Real pres_mbar, amrex::Real &theta, amrex::Real &qv, amrex::Real &qc)
Definition: ERF_SatAdj.H:183
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:274
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:265
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:279
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:275
bool m_do_cond
Definition: ERF_SatAdj.H:276
int m_real_width
Definition: ERF_SatAdj.H:268
@ tabs
Definition: ERF_Kessler.H:24
@ qv
Definition: ERF_Kessler.H:28
@ pres
Definition: ERF_SatAdj.H:33

Referenced by Advance().

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

◆ Copy_Micro_to_State()

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

Updates conserved and microphysics variables in the provided MultiFabs from the internal MultiFabs that store Microphysics module data.

Parameters
[out]consConserved variables
[out]qmoistqv, qc

Reimplemented from NullMoist.

13 {
14  // Get the temperature, density, theta, qt and qp from input
15  for (MFIter mfi(cons,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
16  const auto& tbx = mfi.tilebox();
17 
18  auto states_arr = cons.array(mfi);
19 
20  auto rho_arr = mic_fab_vars[MicVar_SatAdj::rho]->array(mfi);
21  auto theta_arr = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
22  auto qv_arr = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
23  auto qc_arr = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
24 
25  // get potential total density, temperature, qt, qp
26  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
27  {
28  states_arr(i,j,k,RhoTheta_comp) = rho_arr(i,j,k)*theta_arr(i,j,k);
29  states_arr(i,j,k,RhoQ1_comp) = rho_arr(i,j,k)*qv_arr(i,j,k);
30  states_arr(i,j,k,RhoQ2_comp) = rho_arr(i,j,k)*qc_arr(i,j,k);
31  });
32  }
33 
34  // Fill interior ghost cells and periodic boundaries
35  cons.FillBoundary(m_geom.periodicity());
36 }
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
auto qv_arr
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:210
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
@ rho
Definition: ERF_SatAdj.H:30
@ 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 SatAdj::Copy_State_to_Micro ( const amrex::MultiFab &  cons_in)
overridevirtual

Initializes the Microphysics module.

Parameters
[in]cons_inConserved variables input

Reimplemented from NullMoist.

41 {
42  // Get the temperature, density, theta, qt and qp from input
43  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
44  const auto& tbx = mfi.tilebox();
45 
46  auto states_array = cons_in.array(mfi);
47 
48  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
49  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
50 
51  auto rho_array = mic_fab_vars[MicVar_SatAdj::rho]->array(mfi);
52  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
53  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
54  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
55 
56  // Get pressure, theta, temperature, density, and qt, qp
57  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
58  {
59  rho_array(i,j,k) = states_array(i,j,k,Rho_comp);
60  theta_array(i,j,k) = states_array(i,j,k,RhoTheta_comp)/states_array(i,j,k,Rho_comp);
61  qv_array(i,j,k) = states_array(i,j,k,RhoQ1_comp)/states_array(i,j,k,Rho_comp);
62  qc_array(i,j,k) = states_array(i,j,k,RhoQ2_comp)/states_array(i,j,k,Rho_comp);
63 
64  tabs_array(i,j,k) = getTgivenRandRTh(states_array(i,j,k,Rho_comp),
65  states_array(i,j,k,RhoTheta_comp),
66  qv_array(i,j,k));
67 
68  // Pressure in [mbar] for qsat evaluation
69  pres_array(i,j,k) = getPgivenRTh(states_array(i,j,k,RhoTheta_comp), qv_array(i,j,k)) * Real(0.01);
70  });
71  }
72 }
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

Referenced by Update_Micro_Vars().

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

◆ Define()

void SatAdj::Define ( SolverChoice sc)
inlineoverridevirtual

Reimplemented from NullMoist.

58  {
59  m_fac_cond = lcond / sc.c_p;
60  m_rdOcp = sc.rdOcp;
61  m_do_cond = (!sc.use_shoc);
62  }
constexpr amrex::Real lcond
Definition: ERF_Constants.H:77
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1169
amrex::Real c_p
Definition: ERF_DataStruct.H:1168
bool use_shoc
Definition: ERF_DataStruct.H:1200

◆ Init()

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

Initializes the Microphysics module.

Parameters
[in]cons_inConserved variables input
[in]qc_inCloud variables input
[in,out]qv_inVapor variables input
[in]qi_inIce variables input
[in]gridsThe boxes on which we will evolve the solution
[in]geomGeometry associated with these MultiFabs and grids
[in]dt_advanceTimestep for the advance

Reimplemented from NullMoist.

23 {
24  dt = dt_advance;
25  m_geom = geom;
26 
27  // initialize microphysics variables
28  for (auto ivar = 0; ivar < MicVar_SatAdj::NumVars; ++ivar) {
29  mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
30  1, cons_in.nGrowVect());
31  mic_fab_vars[ivar]->setVal(0.);
32  }
33 }
@ NumVars
Definition: ERF_SatAdj.H:37

◆ NewtonSolveSatTemperature()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real SatAdj::NewtonSolveSatTemperature ( const amrex::Real  fac_cond,
const amrex::Real  tabs_old,
const amrex::Real  pres_mbar,
const amrex::Real  qv_old,
amrex::Real qsat 
)
inlinestatic
139  {
140  constexpr amrex::Real tol = amrex::Real(1.0e-8);
141  constexpr int max_niter = 20;
142 
143  amrex::Real dqsat;
144 
145  int niter = 0;
146  amrex::Real fff, dfff;
147  amrex::Real dtabs = one;
148 
149  amrex::Real tabs = tabs_old;
150 
151  //==================================================
152  // Newton iteration to qv=qsat (cloud phase only)
153  //==================================================
154  do {
155  // Saturation moisture fractions
156  erf_qsatw(tabs, pres_mbar, qsat);
157  erf_dtqsatw(tabs, pres_mbar, dqsat);
158 
159  // Function for root finding:
160  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
161  fff = -tabs + tabs_old + fac_cond*(qv_old - qsat);
162 
163  // Derivative of function (T_new iterated on)
164  dfff = -one - fac_cond*dqsat;
165 
166  // Update the temperature
167  dtabs = -fff/dfff;
168  tabs += dtabs;
169 
170  // Update iteration
171  niter = niter+1;
172  } while(std::abs(dtabs) > tol && niter < max_niter);
173 
174  // Update qsat from last iteration (dq = dq/dt * dt)
175  qsat += dqsat*dtabs;
176 
177  return tabs;
178  }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:244

Referenced by AdjustSatAdjCell().

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

◆ Qmoist_Ptr()

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

Reimplemented from NullMoist.

108  {
110  return nullptr;
111  }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
int m_qmoist_size
Definition: ERF_SatAdj.H:259
Here is the call graph for this function:

◆ Qmoist_Restart_Vars()

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

Reimplemented from NullMoist.

123  {
124  a_idx.clear();
125  a_names.clear();
126  }

◆ Qmoist_Size()

int SatAdj::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

114 { return SatAdj::m_qmoist_size; }

◆ Qstate_Moist_Size()

int SatAdj::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

117 { return SatAdj::n_qstate_moist_size; }
int n_qstate_moist_size
Definition: ERF_SatAdj.H:262

◆ Set_RealWidth()

void SatAdj::Set_RealWidth ( const int  real_width)
inlineoverridevirtual

Reimplemented from NullMoist.

129 { m_real_width = real_width; }

◆ Update_Micro_Vars()

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

Reimplemented from NullMoist.

84  {
85  this->Copy_State_to_Micro(cons_in);
86  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:40
Here is the call graph for this function:

◆ Update_State_Vars()

void SatAdj::Update_State_Vars ( amrex::MultiFab &  cons_in,
const amrex::MultiFab &   
)
inlineoverridevirtual

Reimplemented from NullMoist.

92  {
93  this->Copy_Micro_to_State(cons_in);
94  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:12
Here is the call graph for this function:

Member Data Documentation

◆ dt

amrex::Real SatAdj::dt
private

Referenced by Advance().

◆ m_do_cond

bool SatAdj::m_do_cond
private

Referenced by Define().

◆ m_fac_cond

amrex::Real SatAdj::m_fac_cond
private

Referenced by Define().

◆ m_geom

amrex::Geometry SatAdj::m_geom
private

◆ m_qmoist_size

int SatAdj::m_qmoist_size = 0
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_rdOcp

amrex::Real SatAdj::m_rdOcp
private

Referenced by Define().

◆ m_real_width

int SatAdj::m_real_width {0}
private

Referenced by Set_RealWidth().

◆ mic_fab_vars

amrex::Array<FabPtr, MicVar_SatAdj::NumVars> SatAdj::mic_fab_vars
private

◆ n_qstate_moist_size

int SatAdj::n_qstate_moist_size = 2
private

Referenced by Qstate_Moist_Size().


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