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
 
int Qstate_Moist_NumConc_Size () 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
 
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)
 
virtual void Set_RealWidth (const int)
 

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
 
int n_qstate_moist_numconc_size = 0
 
amrex::Geometry m_geom
 
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
51 {}

◆ ~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
211  {
212  amrex::Real qsat;
213  erf_qsatw(tabs, pres_mbar, qsat);
214 
215  if ((qv + qc) > qsat) {
216  // Total water exceeds saturation at the initial temperature, so the final state
217  // is cloudy and saturated. Solve directly for the final saturated temperature.
218  // This is algebraically equivalent to evaporating existing cloud water first
219  // and then recondensing to the fixed-pressure equilibrium state.
220 #if defined(AMREX_DEBUG)
221  const amrex::Real qvprev = qv;
222  const amrex::Real qcprev = qc;
223 #endif
224 
225  if (qc < 0) {
226  // Repair negative cloud water by transferring the deficit back to vapor
227  // before solving the saturated equilibrium and clamping qc to zero.
228  qv += qc;
229  qc = zero;
230  }
231 
232  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
233 
234  const amrex::Real delta_qv = qv - qsat;
235  qv = qsat;
236  qc += delta_qv;
237 
238 #if defined(AMREX_DEBUG)
239  amrex::Real qsatnew;
240  erf_qsatw(tabs, pres_mbar, qsatnew);
241  AMREX_ASSERT(std::abs(qv-qsatnew) < 1e-12);
242  AMREX_ASSERT(std::abs(qv+qc-qvprev-qcprev) < 1e-14);
243 #endif
244 
245  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
246  } else {
247  // Total water does not exceed initial saturation. Evaporate all cloud water,
248  // cooling the cell by latent heat absorption. The cooled cell can become
249  // supersaturated, so recondense if needed.
250  const amrex::Real delta_qc = qc;
251 
252  qv += qc;
253  qc = zero;
254 
255  tabs -= fac_cond * delta_qc;
256  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
257 
258  erf_qsatw(tabs, pres_mbar, qsat);
259  if (qv > qsat) {
260 #if defined(AMREX_DEBUG)
261  const amrex::Real qvprev = qv;
262  const amrex::Real qcprev = qc;
263  const amrex::Real Tprev = tabs;
264 #endif
265 
266  tabs = NewtonSolveSatTemperature(fac_cond, tabs, pres_mbar, qv, qsat);
267 
268  const amrex::Real delta_qv = qv - qsat;
269  qv = qsat;
270  qc += delta_qv;
271 
272 #if defined(AMREX_DEBUG)
273  amrex::Real qsatnew;
274  erf_qsatw(tabs, pres_mbar, qsatnew);
275  AMREX_ASSERT(qv < qvprev);
276  AMREX_ASSERT(qc > qcprev);
277  AMREX_ASSERT(tabs > Tprev);
278  AMREX_ASSERT(std::abs(qv-qsatnew) < 1e-14);
279  AMREX_ASSERT(std::abs(qv+qc-qvprev-qcprev) < 1e-14);
280 #endif
281 
282  theta = getThgivenTandP(tabs, amrex::Real(100.0)*pres_mbar, rdOcp);
283  }
284  }
285  }
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:150
@ theta
Definition: ERF_SatAdj.H:35
@ qv
Definition: ERF_SatAdj.H:39
@ tabs
Definition: ERF_SatAdj.H:36
@ qc
Definition: ERF_SatAdj.H:40
Here is the call graph for this function:

◆ Advance()

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

Reimplemented from NullMoist.

104  {
105  dt = dt_advance;
106 
107  this->AdvanceSatAdj(solverChoice);
108  }
void AdvanceSatAdj(const SolverChoice &)
Definition: ERF_SatAdj.cpp:11
amrex::Real dt
Definition: ERF_SatAdj.H:301
Here is the call graph for this function:

◆ AdvanceSatAdj()

void SatAdj::AdvanceSatAdj ( const SolverChoice )

AdvanceSatAdj is a local valid-cell update over MFIter tileboxes. It uses pressure diagnosed in Copy_State_to_Micro and holds that pressure fixed inside each cell adjustment. There are no stencils, face fluxes, or ghost-cell reads in this kernel.

12 {
13  // Saturation adjustment can be disabled by solver choice, e.g. when SHOC
14  // owns moist thermodynamics instead of the standalone SatAdj module.
15  if (!m_do_cond) { return; }
16 
18 
19  // Expose for GPU
20  Real d_fac_cond = m_fac_cond;
21  Real rdOcp = m_rdOcp;
22 
23  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
24 
25  auto tbx = mfi.tilebox();
26 
27  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
28  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
29  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
30  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
31  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
32 
33  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
34  {
35  Real T = tabs_array(i,j,k);
36  Real p = pres_array(i,j,k);
37  Real th = theta_array(i,j,k);
38  Real qv = qv_array(i,j,k);
39  Real qc = qc_array(i,j,k);
40 
41  AdjustSatAdjCell(d_fac_cond, rdOcp, T, p, th, qv, qc);
42 
43  tabs_array(i,j,k) = T;
44  theta_array(i,j,k) = th;
45  qv_array(i,j,k) = qv;
46  qc_array(i,j,k) = qc;
47  });
48  }
49 }
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
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:204
amrex::Real m_fac_cond
Definition: ERF_SatAdj.H:304
amrex::Array< FabPtr, MicVar_SatAdj::NumVars > mic_fab_vars
Definition: ERF_SatAdj.H:309
amrex::Real m_rdOcp
Definition: ERF_SatAdj.H:305
bool m_do_cond
Definition: ERF_SatAdj.H:306
@ tabs
Definition: ERF_Kessler.H:25
@ qv
Definition: ERF_Kessler.H:29
@ pres
Definition: ERF_SatAdj.H:37
@ p
Definition: ERF_WSM6.H:176

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

Copies SatAdj's adjusted specific variables back into ERF's conserved state. The copied rho is used only to reconstruct rhoTheta, rhoQv, and rhoQc. Density itself is intentionally not written back because SatAdj does not own or modify rho.

Parameters
[out]consConserved variables

Reimplemented from NullMoist.

14 {
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  // Density is intentionally not written back. SatAdj does not modify density;
26  // it uses the copied rho only to form conserved rhoTheta, rhoQv, and rhoQc
27  // from the adjusted specific variables.
28  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
29  {
30  states_arr(i,j,k,RhoTheta_comp) = rho_arr(i,j,k)*theta_arr(i,j,k);
31  states_arr(i,j,k,RhoQ1_comp) = rho_arr(i,j,k)*qv_arr(i,j,k);
32  states_arr(i,j,k,RhoQ2_comp) = rho_arr(i,j,k)*qc_arr(i,j,k);
33  });
34  }
35 
36  // Fill ghost and periodic boundary values after the valid-cell copy-back.
37  cons.FillBoundary(m_geom.periodicity());
38 }
#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
amrex::Geometry m_geom
Definition: ERF_SatAdj.H:298
@ rho
Definition: ERF_SatAdj.H:34
@ cons
Definition: ERF_IndexDefines.H:174

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

Copies conserved ERF state into SatAdj's internal specific variables. rho is copied but not modified by SatAdj. theta, qv, and qc are formed by dividing conserved densities by rho. tabs and pressure are diagnosed from the same conserved state snapshot. pressure is converted from Pa to mbar/hPa for saturation helper calls.

Parameters
[in]cons_inConserved variables input

Reimplemented from NullMoist.

45 {
46  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
47  const auto& tbx = mfi.tilebox();
48 
49  auto states_array = cons_in.array(mfi);
50 
51  auto qv_array = mic_fab_vars[MicVar_SatAdj::qv]->array(mfi);
52  auto qc_array = mic_fab_vars[MicVar_SatAdj::qc]->array(mfi);
53 
54  auto rho_array = mic_fab_vars[MicVar_SatAdj::rho]->array(mfi);
55  auto theta_array = mic_fab_vars[MicVar_SatAdj::theta]->array(mfi);
56  auto tabs_array = mic_fab_vars[MicVar_SatAdj::tabs]->array(mfi);
57  auto pres_array = mic_fab_vars[MicVar_SatAdj::pres]->array(mfi);
58 
59  // Use local scalars so temperature and pressure are diagnosed from the same
60  // conserved state values copied into the microphysics arrays. This avoids
61  // hidden read-after-write dependencies inside the kernel body.
62  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
63  {
64  const Real rho = states_array(i,j,k,Rho_comp);
65  const Real rhoTheta = states_array(i,j,k,RhoTheta_comp);
66  const Real qv = states_array(i,j,k,RhoQ1_comp) / rho;
67  const Real qc = states_array(i,j,k,RhoQ2_comp) / rho;
68  const Real theta = rhoTheta / rho;
69 
70  rho_array(i,j,k) = rho;
71  theta_array(i,j,k) = theta;
72  qv_array(i,j,k) = qv;
73  qc_array(i,j,k) = qc;
74  tabs_array(i,j,k) = getTgivenRandRTh(rho, rhoTheta, qv);
75 
76  // Pressure is stored in mbar/hPa because erf_qsatw expects that unit.
77  pres_array(i,j,k) = getPgivenRTh(rhoTheta, qv) * Real(0.01);
78  });
79  }
80 }
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
rho
Definition: ERF_InitCustomPert_Bubble.H:106
@ theta
Definition: ERF_MM5.H:20

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.

62  {
63  m_fac_cond = lcond / sc.c_p;
64  m_rdOcp = sc.rdOcp;
65  m_do_cond = (!sc.use_shoc);
66  }
constexpr amrex::Real lcond
Definition: ERF_Constants.H:85
amrex::Real rdOcp
Definition: ERF_DataStruct.H:1222
amrex::Real c_p
Definition: ERF_DataStruct.H:1221
bool use_shoc
Definition: ERF_DataStruct.H:1253

◆ 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:41

◆ 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
155  {
156  constexpr amrex::Real tol = amrex::Real(1.0e-8);
157  constexpr int max_niter = 20;
158 
159  amrex::Real dqsat;
160 
161  int niter = 0;
162  amrex::Real fff, dfff;
163  amrex::Real dtabs = one;
164 
165  amrex::Real tabs = tabs_old;
166 
167  //==================================================
168  // Newton iteration to qv=qsat (cloud phase only)
169  //==================================================
170  do {
171  // Saturation moisture fractions
172  erf_qsatw(tabs, pres_mbar, qsat);
173  erf_dtqsatw(tabs, pres_mbar, dqsat);
174 
175  // Function for root finding:
176  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
177  fff = -tabs + tabs_old + fac_cond*(qv_old - qsat);
178 
179  // Derivative of function (T_new iterated on)
180  dfff = -one - fac_cond*dqsat;
181 
182  // Update the temperature
183  dtabs = -fff/dfff;
184  tabs += dtabs;
185 
186  // Update iteration
187  niter = niter+1;
188  } while(std::abs(dtabs) > tol && niter < max_niter);
189 
190  // Return qsat evaluated at the final temperature, not a tangent extrapolation.
191  erf_qsatw(tabs, pres_mbar, qsat);
192 
193  return tabs;
194  }
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.

112  {
114  return nullptr;
115  }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
int m_qmoist_size
Definition: ERF_SatAdj.H:289
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.

130  {
131  a_idx.clear();
132  a_names.clear();
133  }

◆ Qmoist_Size()

int SatAdj::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

118 { return SatAdj::m_qmoist_size; }

◆ Qstate_Moist_NumConc_Size()

int SatAdj::Qstate_Moist_NumConc_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

int n_qstate_moist_numconc_size
Definition: ERF_SatAdj.H:295

◆ Qstate_Moist_Size()

int SatAdj::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

121 { return SatAdj::n_qstate_moist_size; }
int n_qstate_moist_size
Definition: ERF_SatAdj.H:292

◆ Update_Micro_Vars()

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

Reimplemented from NullMoist.

88  {
89  this->Copy_State_to_Micro(cons_in);
90  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSatAdj.cpp:44
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.

96  {
97  this->Copy_Micro_to_State(cons_in);
98  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSatAdj.cpp:13
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().

◆ mic_fab_vars

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

◆ n_qstate_moist_numconc_size

int SatAdj::n_qstate_moist_numconc_size = 0
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: