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

#include <ERF_SAM.H>

Inheritance diagram for SAM:
Collaboration diagram for SAM:

Public Member Functions

 SAM ()
 
virtual ~SAM ()=default
 
void Cloud (const SolverChoice &sc)
 
void IceFall (const SolverChoice &sc)
 
void Precip (const SolverChoice &sc)
 
void PrecipFall (const SolverChoice &sc)
 
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, const amrex::MultiFab &) override
 
void Advance (const amrex::Real &dt_advance, const SolverChoice &sc) override
 
amrex::MultiFab * Qmoist_Ptr (const int &varIdx) override
 
const amrex::MultiFab * Qmoist_Ptr (const int &varIdx) const
 
void Compute_Coefficients ()
 
SAMCoefficientRow CoefficientRowAt (int k) const
 
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_RealWidth (const int)
 

Static Public Member Functions

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat (int &i, int &j, int &k, const int &SAM_moisture_type, const amrex::Real &fac_cond, const amrex::Real &, const amrex::Real &fac_sub, const amrex::Real &an, const amrex::Real &bn, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array, const amrex::Array4< amrex::Real > &qi_array, const amrex::Array4< amrex::Real > &qn_array, const amrex::Array4< amrex::Real > &qt_array)
 

Private Types

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

Private Attributes

int m_qmoist_size = 3
 
int n_qstate_moist_size = 6
 
int n_qstate_moist_numconc_size = 0
 
amrex::Vector< int > MicVarMap
 
amrex::Geometry m_geom
 
amrex::Real dt
 
int nlev
 
int zlo
 
int zhi
 
int m_axis
 
amrex::Real m_fac_cond
 
amrex::Real m_fac_fus
 
amrex::Real m_fac_sub
 
amrex::Real m_rdOcp
 
bool m_do_cond
 
amrex::Real m_dzmin
 
amrex::MultiFab * m_z_phys_nd
 
amrex::MultiFab * m_detJ_cc
 
amrex::Array< FabPtr, MicVar::NumVarsmic_fab_vars
 
amrex::TableData< amrex::Real, 1 > accrrc
 
amrex::TableData< amrex::Real, 1 > accrsi
 
amrex::TableData< amrex::Real, 1 > accrsc
 
amrex::TableData< amrex::Real, 1 > coefice
 
amrex::TableData< amrex::Real, 1 > evaps1
 
amrex::TableData< amrex::Real, 1 > evaps2
 
amrex::TableData< amrex::Real, 1 > accrgi
 
amrex::TableData< amrex::Real, 1 > accrgc
 
amrex::TableData< amrex::Real, 1 > evapg1
 
amrex::TableData< amrex::Real, 1 > evapg2
 
amrex::TableData< amrex::Real, 1 > evapr1
 
amrex::TableData< amrex::Real, 1 > evapr2
 
amrex::TableData< amrex::Real, 1 > rho1d
 
amrex::TableData< amrex::Real, 1 > pres1d
 
amrex::TableData< amrex::Real, 1 > tabs1d
 
amrex::TableData< amrex::Real, 1 > qt1d
 
amrex::TableData< amrex::Real, 1 > qv1d
 
amrex::TableData< amrex::Real, 1 > qn1d
 

Static Private Attributes

static constexpr amrex::Real CFL_MAX = myhalf
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ SAM()

SAM::SAM ( )
inline
60 {}

◆ ~SAM()

virtual SAM::~SAM ( )
virtualdefault

Member Function Documentation

◆ Advance()

void SAM::Advance ( const amrex::Real dt_advance,
const SolverChoice sc 
)
inlineoverridevirtual

Reimplemented from NullMoist.

131  {
132  dt = dt_advance;
133 
134  this->Cloud(sc);
135  this->IceFall(sc);
136  this->Precip(sc);
137  this->PrecipFall(sc);
138  }
void Precip(const SolverChoice &sc)
Definition: ERF_Precip.cpp:11
void IceFall(const SolverChoice &sc)
Definition: ERF_IceFall.cpp:13
void Cloud(const SolverChoice &sc)
Definition: ERF_CloudSAM.cpp:13
void PrecipFall(const SolverChoice &sc)
Definition: ERF_PrecipFall.cpp:17
amrex::Real dt
Definition: ERF_SAM.H:318
Here is the call graph for this function:

◆ Cloud()

void SAM::Cloud ( const SolverChoice sc)

Split cloud components according to saturation pressures; source theta from latent heat.

14 {
15  if (!m_do_cond) { return; }
16 
17  constexpr Real an = one/(tbgmax-tbgmin);
18  constexpr Real bn = tbgmin*an;
19 
20  Real fac_cond = m_fac_cond;
21  Real fac_sub = m_fac_sub;
22  Real fac_fus = m_fac_fus;
23  Real rdOcp = m_rdOcp;
24 
25  int SAM_moisture_type = 1;
26  if (sam_is_no_ice(sc.moisture_type)) {
27  SAM_moisture_type = 2;
28  }
29 
30  for ( MFIter mfi(*(mic_fab_vars[MicVar::tabs]), TilingIfNotGPU()); mfi.isValid(); ++mfi) {
31  auto qt_array = mic_fab_vars[MicVar::qt]->array(mfi);
32  auto qn_array = mic_fab_vars[MicVar::qn]->array(mfi);
33  auto qv_array = mic_fab_vars[MicVar::qv]->array(mfi);
34  auto qcl_array = mic_fab_vars[MicVar::qcl]->array(mfi);
35  auto qci_array = mic_fab_vars[MicVar::qci]->array(mfi);
36 
37  auto tabs_array = mic_fab_vars[MicVar::tabs]->array(mfi);
38  auto theta_array = mic_fab_vars[MicVar::theta]->array(mfi);
39  auto pres_array = mic_fab_vars[MicVar::pres]->array(mfi);
40 
41  auto tbx = mfi.tilebox();
42 
43  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
44  {
45  // Saturation moisture fractions
46  Real omn;
47  Real qsat;
48  Real qsatw;
49  Real qsati;
50 
51  // Newton iteration vars
52  Real delta_qv, delta_qc, delta_qi;
53 
54  // NOTE: Conversion before iterations is necessary to
55  // convert cloud water to ice or vice versa.
56  // This ensures the omn splitting is enforced
57  // before the Newton iteration, which assumes it is.
58 
59  const SAMCloudPhaseChange phase_change =
60  sam_partition_cloud_phase(SAM_moisture_type,
61  tabs_array(i,j,k),
62  qn_array(i,j,k),
63  qcl_array(i,j,k),
64  qci_array(i,j,k),
65  fac_cond, fac_fus,
66  an, bn);
67  omn = phase_change.omn;
68  delta_qc = phase_change.delta_qc;
69  delta_qi = phase_change.delta_qi;
70  qcl_array(i,j,k) = phase_change.qcl;
71  qci_array(i,j,k) = phase_change.qci;
72  tabs_array(i,j,k) = phase_change.tabs;
73  // Cloud adjustment updates tabs under the held-pressure SAM source
74  // convention, then refreshes theta from the same stored pressure.
75  theta_array(i,j,k) = sam_theta_from_stored_mbar_converted_to_pa(tabs_array(i,j,k),
76  pres_array(i,j,k), rdOcp);
77 
78  // Saturation moisture fractions
79  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatw);
80  erf_qsati(tabs_array(i,j,k), pres_array(i,j,k), qsati);
81  qsat = sam_mixed_qsat(omn, qsatw, qsati);
82 
83  // We have enough total moisture to relax to equilibrium
84  if (qt_array(i,j,k) > qsat) {
85 
86  // Update temperature
87  tabs_array(i,j,k) = NewtonIterSat(i, j, k , SAM_moisture_type ,
88  fac_cond , fac_fus , fac_sub ,
89  an , bn ,
90  tabs_array, pres_array,
91  qv_array , qcl_array , qci_array,
92  qn_array , qt_array);
93 
94  // Update theta
95  theta_array(i,j,k) = sam_theta_from_stored_mbar_converted_to_pa(tabs_array(i,j,k),
96  pres_array(i,j,k), rdOcp);
97 
98  //
99  // We cannot blindly relax to qsat, but we can convert qc/qi -> qv.
100  // The concept here is that if we put all the moisture into qv and modify
101  // the temperature, we can then check if qv > qsat occurs (for final T/P/qv).
102  // If the reduction in T/qsat and increase in qv does trigger the
103  // aforementioned condition, we can do Newton iteration to drive qv = qsat.
104  //
105  } else {
106  // Changes in each component
107  delta_qv = qcl_array(i,j,k) + qci_array(i,j,k);
108  delta_qc = qcl_array(i,j,k);
109  delta_qi = qci_array(i,j,k);
110 
111  // Partition the change in non-precipitating q
112  qv_array(i,j,k) += delta_qv;
113  qcl_array(i,j,k) = zero;
114  qci_array(i,j,k) = zero;
115  qn_array(i,j,k) = zero;
116  qt_array(i,j,k) = qv_array(i,j,k);
117 
118  // Update temperature (endothermic since we evap/sublime)
119  tabs_array(i,j,k) -= fac_cond * delta_qc + fac_sub * delta_qi;
120 
121  // Update theta
122  theta_array(i,j,k) = sam_theta_from_stored_mbar_converted_to_pa(tabs_array(i,j,k),
123  pres_array(i,j,k), rdOcp);
124 
125  // Verify assumption that qv > qsat does not occur
126  erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatw);
127  erf_qsati(tabs_array(i,j,k), pres_array(i,j,k), qsati);
128  qsat = sam_mixed_qsat(omn, qsatw, qsati);
129  if (qt_array(i,j,k) > qsat) {
130 
131  // Update temperature
132  tabs_array(i,j,k) = NewtonIterSat(i, j, k , SAM_moisture_type ,
133  fac_cond , fac_fus , fac_sub ,
134  an , bn ,
135  tabs_array, pres_array,
136  qv_array , qcl_array , qci_array,
137  qn_array , qt_array);
138 
139  // Update theta
140  theta_array(i,j,k) = sam_theta_from_stored_mbar_converted_to_pa(tabs_array(i,j,k),
141  pres_array(i,j,k), rdOcp);
142 
143  }
144  }
145  });
146  } // mfi
147 }
constexpr amrex::Real tbgmax
Definition: ERF_Constants.H:51
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real tbgmin
Definition: ERF_Constants.H:50
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_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsati(amrex::Real t, amrex::Real p, amrex::Real &qsati)
Definition: ERF_MicrophysicsUtils.H:218
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 AMREX_FORCE_INLINE SAMCloudPhaseChange sam_partition_cloud_phase(const int SAM_moisture_type, const amrex::Real tabs, const amrex::Real qn, const amrex::Real qcl, const amrex::Real qci, const amrex::Real fac_cond, const amrex::Real fac_fus, const amrex::Real an, const amrex::Real bn) noexcept
Definition: ERF_SAMUtils.H:290
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool sam_is_no_ice(const MoistureType moisture_type) noexcept
Definition: ERF_SAMUtils.H:232
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_mixed_qsat(const amrex::Real omn, const amrex::Real qsatw, const amrex::Real qsati) noexcept
Definition: ERF_SAMUtils.H:343
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_theta_from_stored_mbar_converted_to_pa(const amrex::Real tabs, const amrex::Real pres_mbar, const amrex::Real rdOcp) noexcept
Definition: ERF_SAMUtils.H:171
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Array< FabPtr, MicVar::NumVars > mic_fab_vars
Definition: ERF_SAM.H:341
amrex::Real m_rdOcp
Definition: ERF_SAM.H:330
amrex::Real m_fac_fus
Definition: ERF_SAM.H:328
bool m_do_cond
Definition: ERF_SAM.H:331
amrex::Real m_fac_cond
Definition: ERF_SAM.H:327
amrex::Real m_fac_sub
Definition: ERF_SAM.H:329
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat(int &i, int &j, int &k, const int &SAM_moisture_type, const amrex::Real &fac_cond, const amrex::Real &, const amrex::Real &fac_sub, const amrex::Real &an, const amrex::Real &bn, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array, const amrex::Array4< amrex::Real > &qi_array, const amrex::Array4< amrex::Real > &qn_array, const amrex::Array4< amrex::Real > &qt_array)
Definition: ERF_SAM.H:186
@ pres
Definition: ERF_SAM.H:34
@ qci
Definition: ERF_SAM.H:40
@ qv
Definition: ERF_SAM.H:38
@ qt
Definition: ERF_SAM.H:36
@ qn
Definition: ERF_SAM.H:37
@ theta
Definition: ERF_SAM.H:32
@ qcl
Definition: ERF_SAM.H:39
@ tabs
Definition: ERF_SAM.H:33
@ qsatw
Definition: ERF_WSM6.H:236
@ qsati
Definition: ERF_WSM6.H:236
Definition: ERF_SAMUtils.H:22
amrex::Real delta_qi
Definition: ERF_SAMUtils.H:26
amrex::Real qcl
Definition: ERF_SAMUtils.H:28
amrex::Real tabs
Definition: ERF_SAMUtils.H:32
amrex::Real delta_qc
Definition: ERF_SAMUtils.H:25
amrex::Real omn
Definition: ERF_SAMUtils.H:23
amrex::Real qci
Definition: ERF_SAMUtils.H:29
MoistureType moisture_type
Definition: ERF_DataStruct.H:1299

Referenced by Advance().

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

◆ CoefficientRowAt()

SAMCoefficientRow SAM::CoefficientRowAt ( int  k) const
245 {
246  AMREX_ALWAYS_ASSERT(k >= zlo && k <= zhi);
247 
248  SAMCoefficientRow row{};
249  row.accrrc = copy_table_value_to_host(accrrc, zlo, zhi, k);
250  row.accrsi = copy_table_value_to_host(accrsi, zlo, zhi, k);
251  row.accrsc = copy_table_value_to_host(accrsc, zlo, zhi, k);
252  row.coefice = copy_table_value_to_host(coefice, zlo, zhi, k);
253  row.evaps1 = copy_table_value_to_host(evaps1, zlo, zhi, k);
254  row.evaps2 = copy_table_value_to_host(evaps2, zlo, zhi, k);
255  row.accrgi = copy_table_value_to_host(accrgi, zlo, zhi, k);
256  row.accrgc = copy_table_value_to_host(accrgc, zlo, zhi, k);
257  row.evapg1 = copy_table_value_to_host(evapg1, zlo, zhi, k);
258  row.evapg2 = copy_table_value_to_host(evapg2, zlo, zhi, k);
259  row.evapr1 = copy_table_value_to_host(evapr1, zlo, zhi, k);
260  row.evapr2 = copy_table_value_to_host(evapr2, zlo, zhi, k);
261  return row;
262 }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
int zlo
Definition: ERF_SAM.H:321
amrex::TableData< amrex::Real, 1 > evapr2
Definition: ERF_SAM.H:355
amrex::TableData< amrex::Real, 1 > evapg1
Definition: ERF_SAM.H:352
amrex::TableData< amrex::Real, 1 > accrgi
Definition: ERF_SAM.H:350
amrex::TableData< amrex::Real, 1 > evapg2
Definition: ERF_SAM.H:353
amrex::TableData< amrex::Real, 1 > evaps1
Definition: ERF_SAM.H:348
amrex::TableData< amrex::Real, 1 > accrsi
Definition: ERF_SAM.H:345
amrex::TableData< amrex::Real, 1 > accrrc
Definition: ERF_SAM.H:344
amrex::TableData< amrex::Real, 1 > evapr1
Definition: ERF_SAM.H:354
amrex::TableData< amrex::Real, 1 > evaps2
Definition: ERF_SAM.H:349
amrex::TableData< amrex::Real, 1 > accrgc
Definition: ERF_SAM.H:351
amrex::TableData< amrex::Real, 1 > accrsc
Definition: ERF_SAM.H:346
amrex::TableData< amrex::Real, 1 > coefice
Definition: ERF_SAM.H:347
int zhi
Definition: ERF_SAM.H:321
Definition: ERF_SAMUtils.H:98
amrex::Real accrrc
Definition: ERF_SAMUtils.H:99
Here is the call graph for this function:

◆ Compute_Coefficients()

void SAM::Compute_Coefficients ( )
151 {
152  auto accrrc_t = accrrc.table();
153  auto accrsi_t = accrsi.table();
154  auto accrsc_t = accrsc.table();
155  auto coefice_t = coefice.table();
156  auto evaps1_t = evaps1.table();
157  auto evaps2_t = evaps2.table();
158  auto accrgi_t = accrgi.table();
159  auto accrgc_t = accrgc.table();
160  auto evapg1_t = evapg1.table();
161  auto evapg2_t = evapg2.table();
162  auto evapr1_t = evapr1.table();
163  auto evapr2_t = evapr2.table();
164 
165  auto rho1d_t = rho1d.table();
166  auto pres1d_t = pres1d.table();
167  auto tabs1d_t = tabs1d.table();
168 
169  Real gam3 = erf_gammafff(three );
170  Real gamr1 = erf_gammafff(three+b_rain );
171  Real gamr2 = erf_gammafff((Real(5.0)+b_rain)/two);
172  Real gams1 = erf_gammafff(three+b_snow );
173  Real gams2 = erf_gammafff((Real(5.0)+b_snow)/two);
174  Real gamg1 = erf_gammafff(three+b_grau );
175  Real gamg2 = erf_gammafff((Real(5.0)+b_grau)/two);
176 
177  // calculate the plane average variables
181  rho_ave.compute_averages(ZDir(), rho_ave.field());
182  theta_ave.compute_averages(ZDir(), theta_ave.field());
183  qv_ave.compute_averages(ZDir(), qv_ave.field());
184 
185  // get host variable rho, and rhotheta
186  int ncell = rho_ave.ncell_line();
187 
188  Gpu::HostVector<Real> rho_h(ncell), theta_h(ncell), qv_h(ncell);
189  rho_ave.line_average(0, rho_h);
190  theta_ave.line_average(0, theta_h);
191  qv_ave.line_average(0, qv_h);
192 
193  // copy data to device
194  Gpu::DeviceVector<Real> rho_d(ncell), theta_d(ncell), qv_d(ncell);
195  Gpu::copyAsync(Gpu::hostToDevice, rho_h.begin(), rho_h.end(), rho_d.begin());
196  Gpu::copyAsync(Gpu::hostToDevice, theta_h.begin(), theta_h.end(), theta_d.begin());
197  Gpu::copyAsync(Gpu::hostToDevice, qv_h.begin(), qv_h.end(), qv_d.begin());
198  Gpu::streamSynchronize();
199 
200  Real* rho_dptr = rho_d.data();
201  Real* theta_dptr = theta_d.data();
202  Real* qv_dptr = qv_d.data();
203 
204  ParallelFor(nlev, [=] AMREX_GPU_DEVICE (int k) noexcept
205  {
206  Real RhoTheta = rho_dptr[k]*theta_dptr[k];
207  Real pressure = getPgivenRTh(RhoTheta, qv_dptr[k]);
208  rho1d_t(k) = rho_dptr[k];
209  pres1d_t(k) = sam_pa_to_mbar(pressure);
210  // NOTE: Limit the temperature to the melting point of ice to avoid a divide by
211  // 0 condition when computing the cold evaporation coefficients. This should
212  // not affect results since evaporation requires snow/graupel to be present
213  // and thus T<Real(273.16)
214  tabs1d_t(k) = std::min(getTgivenRandRTh(rho_dptr[k], RhoTheta, qv_dptr[k]),Real(273.16));
215  });
216 
217  if(round(gam3) != 2) {
218  std::cout << "cannot compute gamma-function in Microphysics::Init" << std::endl;
219  std::exit(-1);
220  }
221 
222  // Populate all the coefficients
223  ParallelFor(nlev, [=] AMREX_GPU_DEVICE (int k) noexcept
224  {
225  const SAMCoefficientRow row =
226  sam_compute_coefficient_row(rho1d_t(k), tabs1d_t(k),
227  gamr1, gamr2, gams1, gams2, gamg1, gamg2);
228  accrrc_t(k) = row.accrrc;
229  accrsi_t(k) = row.accrsi;
230  accrsc_t(k) = row.accrsc;
231  coefice_t(k) = row.coefice;
232  evaps1_t(k) = row.evaps1;
233  evaps2_t(k) = row.evaps2;
234  accrgi_t(k) = row.accrgi;
235  accrgc_t(k) = row.accrgc;
236  evapg1_t(k) = row.evapg1;
237  evapg2_t(k) = row.evapg2;
238  evapr1_t(k) = row.evapr1;
239  evapr2_t(k) = row.evapr2;
240  });
241 }
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real b_grau
Definition: ERF_Constants.H:62
constexpr amrex::Real b_rain
Definition: ERF_Constants.H:58
constexpr amrex::Real b_snow
Definition: ERF_Constants.H:60
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
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
pp get("wavelength", wavelength)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real erf_gammafff(amrex::Real x)
Definition: ERF_MicrophysicsUtils.H:19
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMCoefficientRow sam_compute_coefficient_row(const amrex::Real rho, const amrex::Real tabs, const amrex::Real gamr1, const amrex::Real gamr2, const amrex::Real gams1, const amrex::Real gams2, const amrex::Real gamg1, const amrex::Real gamg2) noexcept
Definition: ERF_SAMUtils.H:711
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_pa_to_mbar(const amrex::Real pres_pa) noexcept
Definition: ERF_SAMUtils.H:165
Vector< Real > rho_h(khi+1, zero)
Gpu::DeviceVector< Real > rho_d(khi+1, zero)
Definition: ERF_PlaneAverage.H:14
amrex::TableData< amrex::Real, 1 > rho1d
Definition: ERF_SAM.H:358
int m_axis
Definition: ERF_SAM.H:324
amrex::TableData< amrex::Real, 1 > pres1d
Definition: ERF_SAM.H:359
int nlev
Definition: ERF_SAM.H:321
amrex::Geometry m_geom
Definition: ERF_SAM.H:315
amrex::TableData< amrex::Real, 1 > tabs1d
Definition: ERF_SAM.H:360
@ rho
Definition: ERF_SAM.H:31
amrex::Real coefice
Definition: ERF_SAMUtils.H:102
amrex::Real evaps2
Definition: ERF_SAMUtils.H:104
amrex::Real accrsi
Definition: ERF_SAMUtils.H:100
amrex::Real accrgc
Definition: ERF_SAMUtils.H:106
amrex::Real accrgi
Definition: ERF_SAMUtils.H:105
amrex::Real accrsc
Definition: ERF_SAMUtils.H:101
amrex::Real evapg1
Definition: ERF_SAMUtils.H:107
amrex::Real evapr1
Definition: ERF_SAMUtils.H:109
amrex::Real evapg2
Definition: ERF_SAMUtils.H:108
amrex::Real evapr2
Definition: ERF_SAMUtils.H:110
amrex::Real evaps1
Definition: ERF_SAMUtils.H:103

Referenced by Update_Micro_Vars().

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

◆ Copy_Micro_to_State()

void SAM::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, qi, qr, qs, qg

Reimplemented from NullMoist.

16 {
17  // Get the temperature, density, theta, qt and qp from input
18  for ( MFIter mfi(cons,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
19  const auto& box3d = mfi.tilebox();
20 
21  auto states_arr = cons.array(mfi);
22 
23  auto rho_arr = mic_fab_vars[MicVar::rho]->array(mfi);
24  auto theta_arr = mic_fab_vars[MicVar::theta]->array(mfi);
25 
26  auto qv_arr = mic_fab_vars[MicVar::qv]->array(mfi);
27  auto qc_arr = mic_fab_vars[MicVar::qcl]->array(mfi);
28  auto qi_arr = mic_fab_vars[MicVar::qci]->array(mfi);
29 
30  auto qpr_arr = mic_fab_vars[MicVar::qpr]->array(mfi);
31  auto qps_arr = mic_fab_vars[MicVar::qps]->array(mfi);
32  auto qpg_arr = mic_fab_vars[MicVar::qpg]->array(mfi);
33 
34  // get potential total density, temperature, qt, qp
35  ParallelFor( box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
36  {
37  SAMPrimitiveCell primitive{};
38  primitive.rho = rho_arr(i,j,k);
39  primitive.theta = theta_arr(i,j,k);
40  primitive.qv = qv_arr(i,j,k);
41  primitive.qcl = qc_arr(i,j,k);
42  primitive.qci = qi_arr(i,j,k);
43  primitive.qpr = qpr_arr(i,j,k);
44  primitive.qps = qps_arr(i,j,k);
45  primitive.qpg = qpg_arr(i,j,k);
46  sam_primitive_to_cons(primitive, states_arr, i, j, k);
47  });
48  }
49 
50  // Fill interior ghost cells and periodic boundaries
51  cons.FillBoundary(m_geom.periodicity());
52 }
auto qv_arr
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:210
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void sam_primitive_to_cons(const SAMPrimitiveCell &primitive, const amrex::Array4< amrex::Real > &states_arr, const int i, const int j, const int k) noexcept
Definition: ERF_SAMUtils.H:213
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
@ qpr
Definition: ERF_SAM.H:43
@ qpg
Definition: ERF_SAM.H:45
@ qps
Definition: ERF_SAM.H:44
@ cons
Definition: ERF_IndexDefines.H:174
Definition: ERF_SAMUtils.H:51
amrex::Real rho
Definition: ERF_SAMUtils.H:52

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 SAM::Copy_State_to_Micro ( const amrex::MultiFab &  cons_in)
overridevirtual

Initializes the Microphysics module.

Parameters
[in]cons_inConserved variables input

Reimplemented from NullMoist.

95 {
96  // Get the temperature, density, theta, qt and qp from input
97  for ( MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
98  const auto& box3d = mfi.growntilebox();
99 
100  auto states_array = cons_in.array(mfi);
101 
102  // Non-precipitating
103  auto qv_array = mic_fab_vars[MicVar::qv]->array(mfi);
104  auto qc_array = mic_fab_vars[MicVar::qcl]->array(mfi);
105  auto qi_array = mic_fab_vars[MicVar::qci]->array(mfi);
106  auto qn_array = mic_fab_vars[MicVar::qn]->array(mfi);
107  auto qt_array = mic_fab_vars[MicVar::qt]->array(mfi);
108 
109  // Precipitating
110  auto qpr_array = mic_fab_vars[MicVar::qpr]->array(mfi);
111  auto qps_array = mic_fab_vars[MicVar::qps]->array(mfi);
112  auto qpg_array = mic_fab_vars[MicVar::qpg]->array(mfi);
113  auto qp_array = mic_fab_vars[MicVar::qp]->array(mfi);
114 
115  auto rho_array = mic_fab_vars[MicVar::rho]->array(mfi);
116  auto theta_array = mic_fab_vars[MicVar::theta]->array(mfi);
117  auto tabs_array = mic_fab_vars[MicVar::tabs]->array(mfi);
118  auto pres_array = mic_fab_vars[MicVar::pres]->array(mfi);
119 
120  // Get pressure, theta, temperature, density, and qt, qp
121  ParallelFor(box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
122  {
123  const SAMPrimitiveCell primitive =
124  sam_cons_to_primitive(states_array(i,j,k,Rho_comp),
125  states_array(i,j,k,RhoTheta_comp),
126  states_array(i,j,k,RhoQ1_comp),
127  states_array(i,j,k,RhoQ2_comp),
128  states_array(i,j,k,RhoQ3_comp),
129  states_array(i,j,k,RhoQ4_comp),
130  states_array(i,j,k,RhoQ5_comp),
131  states_array(i,j,k,RhoQ6_comp));
132  rho_array(i,j,k) = primitive.rho;
133  theta_array(i,j,k) = primitive.theta;
134  qv_array(i,j,k) = primitive.qv;
135  qc_array(i,j,k) = primitive.qcl;
136  qi_array(i,j,k) = primitive.qci;
137  qn_array(i,j,k) = primitive.qn;
138  qt_array(i,j,k) = primitive.qt;
139  qpr_array(i,j,k) = primitive.qpr;
140  qps_array(i,j,k) = primitive.qps;
141  qpg_array(i,j,k) = primitive.qpg;
142  qp_array(i,j,k) = primitive.qp;
143  tabs_array(i,j,k) = primitive.tabs;
144  pres_array(i,j,k) = primitive.pres_mbar;
145  });
146  }
147 }
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMPrimitiveCell sam_cons_to_primitive(const amrex::Real rho, const amrex::Real rho_theta, const amrex::Real rho_qv, const amrex::Real rho_qcl, const amrex::Real rho_qci, const amrex::Real rho_qpr, const amrex::Real rho_qps, const amrex::Real rho_qpg) noexcept
Definition: ERF_SAMUtils.H:182
@ qp
Definition: ERF_SAM.H:42
amrex::Real theta
Definition: ERF_SAMUtils.H:53
amrex::Real qpr
Definition: ERF_SAMUtils.H:61
amrex::Real qn
Definition: ERF_SAMUtils.H:59
amrex::Real qps
Definition: ERF_SAMUtils.H:62
amrex::Real qci
Definition: ERF_SAMUtils.H:58
amrex::Real qv
Definition: ERF_SAMUtils.H:56
amrex::Real tabs
Definition: ERF_SAMUtils.H:54
amrex::Real qpg
Definition: ERF_SAMUtils.H:63
amrex::Real qt
Definition: ERF_SAMUtils.H:60
amrex::Real qcl
Definition: ERF_SAMUtils.H:57
amrex::Real pres_mbar
Definition: ERF_SAMUtils.H:55
amrex::Real qp
Definition: ERF_SAMUtils.H:64

Referenced by Update_Micro_Vars().

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

◆ Define()

void SAM::Define ( SolverChoice sc)
inlineoverridevirtual

Reimplemented from NullMoist.

80  {
81  m_fac_cond = lcond / sc.c_p;
82  m_fac_fus = lfus / sc.c_p;
83  m_fac_sub = lsub / sc.c_p;
84  m_axis = sc.ave_plane;
85  m_rdOcp = sc.rdOcp;
86  m_do_cond = (!sc.use_shoc);
87  }
constexpr amrex::Real lfus
Definition: ERF_Constants.H:86
constexpr amrex::Real lsub
Definition: ERF_Constants.H:87
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
int ave_plane
Definition: ERF_DataStruct.H:1313

◆ IceFall()

void SAM::IceFall ( const SolverChoice sc)

Sedimentation of cloud ice (A32)

13  {
14 
16  return;
17 
18  Real dtn = dt;
19  Real coef = dtn/m_dzmin;
20 
21  auto domain = m_geom.Domain();
22  int k_lo = domain.smallEnd(2);
23  int k_hi = domain.bigEnd(2);
24 
27  auto qn = mic_fab_vars[MicVar::qn];
28  auto qt = mic_fab_vars[MicVar::qt];
31 
32  MultiFab fz;
33  IntVect ng = qcl->nGrowVect();
34  BoxArray ba = qcl->boxArray();
35  DistributionMapping dm = qcl->DistributionMap();
36  fz.define(convert(ba, IntVect(0,0,1)), dm, 1, ng);
37  fz.setVal(0.);
38 
39  for (MFIter mfi(fz, TileNoZ()); mfi.isValid(); ++mfi) {
40  auto qci_array = qci->array(mfi);
41  auto rho_array = rho->array(mfi);
42  auto fz_array = fz.array(mfi);
43 
44  const auto& box3d = mfi.tilebox();
45 
46  ParallelFor(box3d, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
47  {
48  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
49  const Real rho_k = (k == k_hi + 1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
50  const Real qci_km1 = (k == k_lo) ? qci_array(i,j,k) : qci_array(i,j,k-1);
51  const Real qci_k = (k == k_hi + 1) ? qci_array(i,j,k-1) : qci_array(i,j,k);
52  const SAMFaceState face_state = sam_face_average_state(k, k_lo, k_hi,
53  rho_km1, rho_k,
54  zero, zero,
55  qci_km1, qci_k,
56  zero, zero);
57  const Real vt_ice = sam_cloud_ice_terminal_velocity(face_state.qci_avg);
58 
59  // NOTE: Fz is the sedimentation flux from the advective operator.
60  // In the terrain-following coordinate system, the z-deriv in
61  // the divergence uses the normal velocity (Omega). However,
62  // there are no u/v components to the sedimentation velocity.
63  // Therefore, we simply end up with a division by detJ when
64  // evaluating the source term: dJinv * (flux_hi - flux_lo) * dzinv.
65  fz_array(i,j,k) = face_state.rho_avg * vt_ice * face_state.qci_avg;
66  });
67  }
68 
69  // Compute number of substeps from maximum terminal velocity
70  Real wt_max;
71  int n_substep;
72  auto const& ma_fz_arr = fz.const_arrays();
73  GpuTuple<Real> max = ParReduce(TypeList<ReduceOpMax>{},
74  TypeList<Real>{},
75  fz, IntVect(0),
76  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
77  -> GpuTuple<Real>
78  {
79  return { ma_fz_arr[box_no](i,j,k) };
80  });
81  wt_max = get<0>(max);
82  // ParReduce over a MultiFab gives the local rank maximum here. IceFall
83  // needs one global substep count across ranks, matching PrecipFall.
84  ParallelDescriptor::ReduceRealMax(wt_max);
86  n_substep = sam_substep_count_from_reduced_flux(wt_max, dtn, m_dzmin);
87  AMREX_ALWAYS_ASSERT(n_substep >= 1);
88  coef /= Real(n_substep);
89  dtn /= Real(n_substep);
90 
91  // Substep the vertical advection
92  for (int nsub(0); nsub<n_substep; ++nsub) {
93  for (MFIter mfi(*qci, TileNoZ()); mfi.isValid(); ++mfi) {
94  auto qci_array = qci->array(mfi);
95  auto qn_array = qn->array(mfi);
96  auto qt_array = qt->array(mfi);
97  auto rho_array = rho->array(mfi);
98  auto fz_array = fz.array(mfi);
99 
100  const auto dJ_array = (m_detJ_cc) ? m_detJ_cc->const_array(mfi) : Array4<const Real>{};
101 
102  const auto& tbx = mfi.tilebox();
103  const auto& tbz = mfi.tilebox(IntVect(0,0,1),IntVect(0));
104 
105  // Update vertical flux every substep
106  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
107  {
108  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
109  const Real rho_k = (k == k_hi + 1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
110  const Real qci_km1 = (k == k_lo) ? qci_array(i,j,k) : qci_array(i,j,k-1);
111  const Real qci_k = (k == k_hi + 1) ? qci_array(i,j,k-1) : qci_array(i,j,k);
112  const SAMFaceState face_state = sam_face_average_state(k, k_lo, k_hi,
113  rho_km1, rho_k,
114  zero, zero,
115  qci_km1, qci_k,
116  zero, zero);
117  const Real vt_ice = sam_cloud_ice_terminal_velocity(face_state.qci_avg);
118 
119  // NOTE: Fz is the sedimentation flux from the advective operator.
120  // In the terrain-following coordinate system, the z-deriv in
121  // the divergence uses the normal velocity (Omega). However,
122  // there are no u/v components to the sedimentation velocity.
123  // Therefore, we simply end up with a division by detJ when
124  // evaluating the source term: dJinv * (flux_hi - flux_lo) * dzinv.
125  fz_array(i,j,k) = face_state.rho_avg * vt_ice * face_state.qci_avg;
126  });
127 
128  // Update precip every substep
129  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
130  {
131  // Jacobian determinant
132  Real dJinv = (dJ_array) ? one/dJ_array(i,j,k) : one;
133 
134  //==================================================
135  // Cloud ice sedimentation (A32)
136  //==================================================
137  Real dqi = sam_sedimentation_tendency(fz_array(i,j,k+1), fz_array(i,j,k),
138  rho_array(i,j,k), dJinv, coef);
139  dqi = std::max(-qci_array(i,j,k), dqi);
140 
141  // Add this increment to both non-precipitating and total water.
142  qci_array(i,j,k) += dqi;
143  qn_array(i,j,k) += dqi;
144  qt_array(i,j,k) += dqi;
145 
146  // NOTE: Sedimentation does not affect the potential temperature,
147  // but it does affect the liquid/ice static energy.
148  // No source to Theta occurs here.
149  });
150  } // mfi
151  } // nsub
152 }
rho
Definition: ERF_InitCustomPert_Bubble.H:106
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int sam_substep_count_from_reduced_flux(const amrex::Real reduced_flux, const amrex::Real dt, const amrex::Real dz) noexcept
Definition: ERF_SAMUtils.H:823
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_cloud_ice_terminal_velocity(const amrex::Real qci_avg) noexcept
Definition: ERF_SAMUtils.H:768
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_sedimentation_tendency(const amrex::Real fz_hi, const amrex::Real fz_lo, const amrex::Real rho, const amrex::Real dJinv, const amrex::Real coef) noexcept
Definition: ERF_SAMUtils.H:1015
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMFaceState sam_face_average_state(const int k, const int k_lo, const int k_hi, const amrex::Real rho_km1, const amrex::Real rho_k, const amrex::Real tabs_km1, const amrex::Real tabs_k, const amrex::Real qci_km1, const amrex::Real qci_k, const amrex::Real qp_km1, const amrex::Real qp_k) noexcept
Definition: ERF_SAMUtils.H:776
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
amrex::MultiFab * m_detJ_cc
Definition: ERF_SAM.H:338
amrex::Real m_dzmin
Definition: ERF_SAM.H:334
@ qcl
Definition: ERF_Kessler.H:30
@ tabs
Definition: ERF_Kessler.H:25
@ qt
Definition: ERF_Kessler.H:28
@ ng
Definition: ERF_Morrison.H:48
@ qn
Definition: ERF_Morrison.H:33
@ qci
Definition: ERF_Morrison.H:36
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
Definition: ERF_SAMUtils.H:132
amrex::Real qci_avg
Definition: ERF_SAMUtils.H:135
amrex::Real rho_avg
Definition: ERF_SAMUtils.H:133

Referenced by Advance().

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

◆ Init()

void SAM::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

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.

43 {
44  dt = dt_advance;
45  m_geom = geom;
46 
47  m_z_phys_nd = z_phys_nd.get();
48  m_detJ_cc = detJ_cc.get();
49 
50  MicVarMap.resize(m_qmoist_size);
52 
53  // initialize microphysics variables
54  for (auto ivar = 0; ivar < MicVar::NumVars; ++ivar) {
55  mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
56  1, cons_in.nGrowVect());
57  mic_fab_vars[ivar]->setVal(0.);
58  }
59 
60  // NOTE: For multi-level not all ranks will own a box.
61  // Furthermore, the plane average allocates space
62  // for the entire domain. We make this consistent.
63  nlev = m_geom.Domain().length(2);
64  zlo = m_geom.Domain().smallEnd(2);
65  zhi = m_geom.Domain().bigEnd(2);
66 
67  // parameters
68  accrrc.resize({zlo}, {zhi});
69  accrsi.resize({zlo}, {zhi});
70  accrsc.resize({zlo}, {zhi});
71  coefice.resize({zlo}, {zhi});
72  evaps1.resize({zlo}, {zhi});
73  evaps2.resize({zlo}, {zhi});
74  accrgi.resize({zlo}, {zhi});
75  accrgc.resize({zlo}, {zhi});
76  evapg1.resize({zlo}, {zhi});
77  evapg2.resize({zlo}, {zhi});
78  evapr1.resize({zlo}, {zhi});
79  evapr2.resize({zlo}, {zhi});
80 
81  // data (input)
82  rho1d.resize({zlo}, {zhi});
83  pres1d.resize({zlo}, {zhi});
84  tabs1d.resize({zlo}, {zhi});
85 }
int m_qmoist_size
Definition: ERF_SAM.H:300
amrex::MultiFab * m_z_phys_nd
Definition: ERF_SAM.H:337
amrex::Vector< int > MicVarMap
Definition: ERF_SAM.H:312
@ NumVars
Definition: ERF_SAM.H:50
@ snow_accum
Definition: ERF_SAM.H:48
@ rain_accum
Definition: ERF_SAM.H:47
@ graup_accum
Definition: ERF_SAM.H:49

◆ NewtonIterSat()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real SAM::NewtonIterSat ( int &  i,
int &  j,
int &  k,
const int &  SAM_moisture_type,
const amrex::Real fac_cond,
const amrex::Real ,
const amrex::Real fac_sub,
const amrex::Real an,
const amrex::Real bn,
const amrex::Array4< amrex::Real > &  tabs_array,
const amrex::Array4< amrex::Real > &  pres_array,
const amrex::Array4< amrex::Real > &  qv_array,
const amrex::Array4< amrex::Real > &  qc_array,
const amrex::Array4< amrex::Real > &  qi_array,
const amrex::Array4< amrex::Real > &  qn_array,
const amrex::Array4< amrex::Real > &  qt_array 
)
inlinestatic
200  {
201  // Solution tolerance
202  amrex::Real tol = amrex::Real(1.0e-4);
203 
204  // Saturation moisture fractions
205  amrex::Real omn, domn;
206  amrex::Real qsat, dqsat;
207  amrex::Real qsatw, dqsatw;
208  amrex::Real qsati, dqsati;
209 
210  // Newton iteration vars
211  int niter;
212  amrex::Real fff, dfff, dtabs;
213  amrex::Real lstar, dlstar;
214  amrex::Real lstarw, lstari;
215  amrex::Real delta_qc, delta_qi;
216 
217  // Initial guess for temperature & pressure
218  amrex::Real tabs = tabs_array(i,j,k);
219  amrex::Real pres = pres_array(i,j,k);
220 
221  niter = 0;
222  dtabs = 1;
223  //==================================================
224  // Newton iteration to qv=qsat (cloud phase only)
225  //==================================================
226  do {
227  // Latent heats and their derivatives wrt to T
228  lstarw = fac_cond;
229  lstari = fac_sub;
230  domn = zero;
231 
232  // Saturation moisture fractions
235  erf_dtqsatw(tabs, pres, dqsatw);
236  erf_dtqsati(tabs, pres, dqsati);
237 
238  omn = sam_cloud_liquid_fraction(SAM_moisture_type, tabs, an, bn);
239 
240  if (SAM_moisture_type == 1) {
241  // Cloud ice not permitted (condensation)
242  if(tabs >= tbgmax) {
243  }
244  // Cloud water not permitted (sublimation)
245  else if(tabs <= tbgmin) {
246  }
247  // Mixed cloud phase (condensation & deposition)
248  else {
249  lstari = fac_sub;
250  domn = an;
251  }
252  } else if (SAM_moisture_type == 2) {
253  domn = zero;
254  }
255 
256  // Linear combination of each component
257  qsat = sam_mixed_qsat(omn, qsatw, qsati);
258  dqsat = sam_mixed_dqsat_dT(omn, domn, qsatw, qsati, dqsatw, dqsati);
259  lstar = omn * lstarw + (one-omn) * lstari;
260  dlstar = domn * lstarw - domn * lstari;
261 
262  // Function for root finding:
263  // 0 = -T_new + T_old + L_eff/C_p * (qv - qsat)
264  fff = sam_newton_residual(tabs, tabs_array(i,j,k), lstar, qv_array(i,j,k), qsat);
265 
266  // Derivative of function (T_new iterated on)
267  dfff = sam_newton_residual_derivative(lstar, dlstar, qv_array(i,j,k), qsat, dqsat);
268 
269  // Update the temperature
270  dtabs = -fff/dfff;
271  tabs += dtabs;
272 
273  // Update iteration
274  niter = niter+1;
275  } while(std::abs(dtabs) > tol && niter < 20);
276 
277  // Update qsat from last iteration (dq = dq/dt * dt)
278  qsat += dqsat*dtabs;
279 
280  // Changes in each component
281  const SAMCloudPhaseChange phase_change =
282  sam_apply_condensate_limiter(qv_array(i,j,k), qsat,
283  qc_array(i,j,k), qi_array(i,j,k), omn);
284  delta_qc = phase_change.delta_qc;
285  delta_qi = phase_change.delta_qi;
286 
287  // Partition the change in non-precipitating q
288  qv_array(i,j,k) = phase_change.qv;
289  qc_array(i,j,k) += delta_qc;
290  qi_array(i,j,k) += delta_qi;
291  qn_array(i,j,k) = qc_array(i,j,k) + qi_array(i,j,k);
292  qt_array(i,j,k) = qv_array(i,j,k) + qn_array(i,j,k);
293 
294  // Return to temperature
295  return tabs;
296  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:244
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsati(amrex::Real t, amrex::Real p, amrex::Real &dtqsati)
Definition: ERF_MicrophysicsUtils.H:235
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_newton_residual_derivative(const amrex::Real lstar, const amrex::Real dlstar, const amrex::Real qv, const amrex::Real qsat, const amrex::Real dqsat) noexcept
Definition: ERF_SAMUtils.H:376
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_cloud_liquid_fraction(const int SAM_moisture_type, const amrex::Real tabs, const amrex::Real an, const amrex::Real bn) noexcept
Definition: ERF_SAMUtils.H:245
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_mixed_dqsat_dT(const amrex::Real omn, const amrex::Real domn, const amrex::Real qsatw, const amrex::Real qsati, const amrex::Real dqsatw, const amrex::Real dqsati) noexcept
Definition: ERF_SAMUtils.H:352
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMCloudPhaseChange sam_apply_condensate_limiter(const amrex::Real qv, const amrex::Real qsat, const amrex::Real qc, const amrex::Real qi, const amrex::Real omn) noexcept
Definition: ERF_SAMUtils.H:387
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_newton_residual(const amrex::Real tabs_new, const amrex::Real tabs_old, const amrex::Real lstar, const amrex::Real qv, const amrex::Real qsat) noexcept
Definition: ERF_SAMUtils.H:365
amrex::Real qv
Definition: ERF_SAMUtils.H:27
Here is the call graph for this function:

◆ Precip()

void SAM::Precip ( const SolverChoice sc)

Autoconversion (A30), Accretion (A28), Evaporation (A24)

12 {
13  if (sam_is_no_precip(sc.moisture_type)) return;
14 
15  Real powr1 = (three + b_rain) / Real(4.0);
16  Real powr2 = (Real(5.0) + b_rain) / Real(8.0);
17  Real pows1 = (three + b_snow) / Real(4.0);
18  Real pows2 = (Real(5.0) + b_snow) / Real(8.0);
19  Real powg1 = (three + b_grau) / Real(4.0);
20  Real powg2 = (Real(5.0) + b_grau) / Real(8.0);
21 
22  auto accrrc_t = accrrc.table();
23  auto accrsc_t = accrsc.table();
24  auto accrsi_t = accrsi.table();
25  auto accrgc_t = accrgc.table();
26  auto accrgi_t = accrgi.table();
27  auto coefice_t = coefice.table();
28  auto evapr1_t = evapr1.table();
29  auto evapr2_t = evapr2.table();
30  auto evaps1_t = evaps1.table();
31  auto evaps2_t = evaps2.table();
32  auto evapg1_t = evapg1.table();
33  auto evapg2_t = evapg2.table();
34 
35  Real fac_cond = m_fac_cond;
36  Real fac_sub = m_fac_sub;
37  Real fac_fus = m_fac_fus;
38  Real rdOcp = m_rdOcp;
39 
41 
42  Real dtn = dt;
43 
44  int SAM_moisture_type = 1;
45  if (sc.moisture_type == MoistureType::SAM_NoIce) {
46  SAM_moisture_type = 2;
47  }
48 
49  const SAMPrecipConfig precip_config{
50  SAM_moisture_type,
51  true,
52  dtn,
53  rdOcp,
54  fac_cond,
55  fac_fus,
56  fac_sub,
57  eps,
58  powr1,
59  pows1,
60  powg1,
61  powr2,
62  pows2,
63  powg2};
64  // get the temperature, density, theta, qt and qp from input
65  for ( MFIter mfi(*(mic_fab_vars[MicVar::tabs]),TilingIfNotGPU()); mfi.isValid(); ++mfi) {
66  auto rho_array = mic_fab_vars[MicVar::rho]->array(mfi);
67  auto theta_array = mic_fab_vars[MicVar::theta]->array(mfi);
68  auto tabs_array = mic_fab_vars[MicVar::tabs]->array(mfi);
69  auto pres_array = mic_fab_vars[MicVar::pres]->array(mfi);
70 
71  // Non-precipitating
72  auto qv_array = mic_fab_vars[MicVar::qv]->array(mfi);
73  auto qcl_array = mic_fab_vars[MicVar::qcl]->array(mfi);
74  auto qci_array = mic_fab_vars[MicVar::qci]->array(mfi);
75  auto qn_array = mic_fab_vars[MicVar::qn]->array(mfi);
76  auto qt_array = mic_fab_vars[MicVar::qt]->array(mfi);
77 
78  // Precipitating
79  auto qpr_array = mic_fab_vars[MicVar::qpr]->array(mfi);
80  auto qps_array = mic_fab_vars[MicVar::qps]->array(mfi);
81  auto qpg_array = mic_fab_vars[MicVar::qpg]->array(mfi);
82  auto qp_array = mic_fab_vars[MicVar::qp]->array(mfi);
83 
84  auto tbx = mfi.tilebox();
85 
86  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
87  {
88  SAMCellState state{};
89  state.rho = rho_array(i,j,k);
90  state.theta = theta_array(i,j,k);
91  state.tabs = tabs_array(i,j,k);
92  state.pres_mbar = pres_array(i,j,k);
93  state.qv = qv_array(i,j,k);
94  state.qcl = qcl_array(i,j,k);
95  state.qci = qci_array(i,j,k);
96  state.qn = qn_array(i,j,k);
97  state.qt = qt_array(i,j,k);
98  state.qpr = qpr_array(i,j,k);
99  state.qps = qps_array(i,j,k);
100  state.qpg = qpg_array(i,j,k);
101  state.qp = qp_array(i,j,k);
102 
103  const SAMCoefficientRow coeffs{
104  accrrc_t(k),
105  accrsi_t(k),
106  accrsc_t(k),
107  coefice_t(k),
108  evaps1_t(k),
109  evaps2_t(k),
110  accrgi_t(k),
111  accrgc_t(k),
112  evapg1_t(k),
113  evapg2_t(k),
114  evapr1_t(k),
115  evapr2_t(k)};
116 
117  // The scalar helper owns the local held-pressure source update;
118  // this public kernel handles state staging and coefficient lookup.
119  state = sam_precip_cell_update(state, coeffs, precip_config);
120 
121  theta_array(i,j,k) = state.theta;
122  tabs_array(i,j,k) = state.tabs;
123  qv_array(i,j,k) = state.qv;
124  qcl_array(i,j,k) = state.qcl;
125  qci_array(i,j,k) = state.qci;
126  qn_array(i,j,k) = state.qn;
127  qt_array(i,j,k) = state.qt;
128  qpr_array(i,j,k) = state.qpr;
129  qps_array(i,j,k) = state.qps;
130  qpg_array(i,j,k) = state.qpg;
131  qp_array(i,j,k) = state.qp;
132  });
133  }
134 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool sam_is_no_precip(const MoistureType moisture_type) noexcept
Definition: ERF_SAMUtils.H:239
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMCellState sam_precip_cell_update(SAMCellState state, const SAMCoefficientRow &coeffs, const SAMPrecipConfig &config, SAMPrecipCellDiagnostics *diagnostics=nullptr) noexcept
Definition: ERF_SAMUtils.H:559
Definition: ERF_SAMUtils.H:69

Referenced by Advance().

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

◆ PrecipFall()

void SAM::PrecipFall ( const SolverChoice sc)

Precipitation fluxes P_{r/s/g} (A19)

Code modified from SAMXX, the C++ version of the SAM code.

Parameters
[in]hydro_typeType selection for the precipitation advection hydrodynamics scheme (0-3)
18 {
19  if (sam_is_no_precip(sc.moisture_type)) return;
20 
21  // Conservative component-sedimentation contract: rain, snow, and graupel
22  // sediment with separate component fluxes. Each face flux is limited by
23  // the donor cell's available component mass, including detJ when present.
24  // The same limited bottom fluxes update surface accumulation and the
25  // in-column qpr/qps/qpg flux-divergence update. This closes component and
26  // total precipitation column budgets up to roundoff. Sedimentation does
27  // not update theta.
28 
29  Real rho_0 = Real(1.29);
30 
31  Real gamr3 = erf_gammafff(Real(4.0)+b_rain);
32  Real gams3 = erf_gammafff(Real(4.0)+b_snow);
33  Real gamg3 = erf_gammafff(Real(4.0)+b_grau);
34 
35  Real vrain = (a_rain*gamr3/Real(6.0))*std::pow((PI*rhor*nzeror),-crain);
36  Real vsnow = (a_snow*gams3/Real(6.0))*std::pow((PI*rhos*nzeros),-csnow);
37  Real vgrau = (a_grau*gamg3/Real(6.0))*std::pow((PI*rhog*nzerog),-cgrau);
38 
39  Real dtn = dt;
40  Real coef = dtn/m_dzmin;
41 
42  auto domain = m_geom.Domain();
43  int k_lo = domain.smallEnd(2);
44  int k_hi = domain.bigEnd(2);
45 
49  auto qp = mic_fab_vars[MicVar::qp];
55 
56  auto ba = tabs->boxArray();
57  auto dm = tabs->DistributionMap();
58  auto ngrow = tabs->nGrowVect();
59 
60  MultiFab fz;
61  fz.define(convert(ba, IntVect(0,0,1)), dm, 3, ngrow);
62 
63  // Precompute the vertical component fluxes for the legacy reduced-flux substep rule.
64  for (MFIter mfi(fz, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
65  auto qpr_array = qpr->const_array(mfi);
66  auto qps_array = qps->const_array(mfi);
67  auto qpg_array = qpg->const_array(mfi);
68  auto rho_array = rho->const_array(mfi);
69  auto tabs_array = tabs->const_array(mfi);
70  auto fz_array = fz.array(mfi);
71 
72  const auto& box3d = mfi.tilebox();
73 
74  ParallelFor(box3d, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
75  {
76  const SAMPrecipComponentFaceState face_state =
77  sam_precip_component_face_state(rho_array, tabs_array,
78  qpr_array, qps_array, qpg_array,
79  i, j, k, k_lo, k_hi);
80  const SAMPrecipFluxComponents raw_fluxes =
82  vrain, vsnow, vgrau);
83  const SAMPrecipFluxComponents corrected_fluxes =
85  rho_0,
86  face_state.rho_avg);
87 
88  fz_array(i,j,k,0) = corrected_fluxes.rain;
89  fz_array(i,j,k,1) = corrected_fluxes.snow;
90  fz_array(i,j,k,2) = corrected_fluxes.graupel;
91  });
92  }
93 
94  // Compute the legacy reduced-flux substep count from the maximum
95  // density-corrected sedimentation flux rather than a direct fall speed.
96  Real max_reduced_flux;
97  int n_substep;
98  auto const& ma_fz_arr = fz.const_arrays();
99  GpuTuple<Real> max = ParReduce(TypeList<ReduceOpMax>{},
100  TypeList<Real>{},
101  fz, IntVect(0),
102  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
103  -> GpuTuple<Real>
104  {
105  return { ma_fz_arr[box_no](i,j,k,0)
106  + ma_fz_arr[box_no](i,j,k,1)
107  + ma_fz_arr[box_no](i,j,k,2) };
108  });
109  max_reduced_flux = get<0>(max);
110  // ParReduce over a MultiFab gives the local rank maximum here. PrecipFall
111  // needs one global substep count across ranks, including ranks that own no
112  // active tiles for this reduction, so keep the explicit MPI max reduction.
113  ParallelDescriptor::ReduceRealMax(max_reduced_flux);
114  max_reduced_flux += std::numeric_limits<Real>::epsilon();
115  n_substep = sam_substep_count_from_reduced_flux(max_reduced_flux, dtn, m_dzmin);
116  AMREX_ALWAYS_ASSERT(n_substep >= 1);
117  coef /= Real(n_substep);
118  dtn /= Real(n_substep);
119 
120  // Substep the vertical advection
121  for (int nsub(0); nsub<n_substep; ++nsub) {
122  for (MFIter mfi(*qp, TileNoZ()); mfi.isValid(); ++mfi) {
123  auto qpr_array = qpr->array(mfi);
124  auto qps_array = qps->array(mfi);
125  auto qpg_array = qpg->array(mfi);
126  auto qp_array = qp->array(mfi);
127  auto rho_array = rho->const_array(mfi);
128  auto tabs_array = tabs->const_array(mfi);
129  auto fz_array = fz.array(mfi);
130 
131  auto rain_accum_array = rain_accum->array(mfi);
132  auto snow_accum_array = snow_accum->array(mfi);
133  auto graup_accum_array = graup_accum->array(mfi);
134 
135  const auto dJ_array = (m_detJ_cc) ? m_detJ_cc->const_array(mfi) : Array4<const Real>{};
136 
137  const auto& tbx = mfi.tilebox();
138  const auto& tbz = mfi.tilebox(IntVect(0,0,1),IntVect(0));
139 
140  // Update vertical flux every substep
141  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
142  {
143  const SAMPrecipComponentFaceState face_state =
144  sam_precip_component_face_state(rho_array, tabs_array,
145  qpr_array, qps_array, qpg_array,
146  i, j, k, k_lo, k_hi);
147  const SAMPrecipFluxComponents raw_fluxes =
149  vrain, vsnow, vgrau);
150  const SAMPrecipFluxComponents corrected_fluxes =
152  rho_0,
153  face_state.rho_avg);
154  const int donor_k = sam_precip_face_donor_k(k, k_lo, k_hi);
155  const Real detJ_donor = (dJ_array) ? dJ_array(i,j,donor_k) : one;
156 
157  // NOTE: Fz is the sedimentation flux from the advective operator.
158  // In the terrain-following coordinate system, the z-deriv in
159  // the divergence uses the normal velocity (Omega). However,
160  // there are no u/v components to the sedimentation velocity.
161  // Therefore, we simply end up with a division by detJ when
162  // evaluating the source term: dJinv * (flux_hi - flux_lo) * dzinv.
163  fz_array(i,j,k,0) = sam_limit_precip_component_flux(corrected_fluxes.rain,
164  rho_array(i,j,donor_k),
165  qpr_array(i,j,donor_k),
166  detJ_donor,
167  coef);
168  fz_array(i,j,k,1) = sam_limit_precip_component_flux(corrected_fluxes.snow,
169  rho_array(i,j,donor_k),
170  qps_array(i,j,donor_k),
171  detJ_donor,
172  coef);
173  fz_array(i,j,k,2) = sam_limit_precip_component_flux(corrected_fluxes.graupel,
174  rho_array(i,j,donor_k),
175  qpg_array(i,j,donor_k),
176  detJ_donor,
177  coef);
178 
179  if (k == k_lo) {
180  const SAMSurfaceAccumulation surface_accum =
182  {fz_array(i,j,k,0), fz_array(i,j,k,1), fz_array(i,j,k,2)},
183  dtn);
184  rain_accum_array(i,j,k) = rain_accum_array(i,j,k) + surface_accum.rain;
185  snow_accum_array(i,j,k) = snow_accum_array(i,j,k) + surface_accum.snow;
186  graup_accum_array(i,j,k) = graup_accum_array(i,j,k) + surface_accum.graupel;
187  }
188  });
189 
190  // Update precip every substep
191  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
192  {
193  // Jacobian determinant
194  Real dJinv = (dJ_array) ? one/dJ_array(i,j,k) : one;
195 
196  //==================================================
197  // Precipitating sedimentation (A19)
198  //==================================================
199  const Real dqpr = sam_sedimentation_tendency(fz_array(i,j,k+1,0), fz_array(i,j,k,0),
200  rho_array(i,j,k), dJinv, coef);
201  const Real dqps = sam_sedimentation_tendency(fz_array(i,j,k+1,1), fz_array(i,j,k,1),
202  rho_array(i,j,k), dJinv, coef);
203  const Real dqpg = sam_sedimentation_tendency(fz_array(i,j,k+1,2), fz_array(i,j,k,2),
204  rho_array(i,j,k), dJinv, coef);
205 
206  qpr_array(i,j,k) = std::max(Real(0.0), qpr_array(i,j,k) + dqpr);
207  qps_array(i,j,k) = std::max(Real(0.0), qps_array(i,j,k) + dqps);
208  qpg_array(i,j,k) = std::max(Real(0.0), qpg_array(i,j,k) + dqpg);
209  qp_array(i,j,k) = qpr_array(i,j,k) + qps_array(i,j,k) + qpg_array(i,j,k);
210 
211  // NOTE: Sedimentation does not affect the potential temperature,
212  // but it does affect the liquid/ice static energy.
213  // No source to Theta occurs here.
214  });
215  } // mfi
216  } // nsub
217 }
constexpr amrex::Real csnow
Definition: ERF_Constants.H:101
constexpr amrex::Real rhog
Definition: ERF_Constants.H:49
constexpr amrex::Real nzerog
Definition: ERF_Constants.H:78
constexpr amrex::Real cgrau
Definition: ERF_Constants.H:102
constexpr amrex::Real a_grau
Definition: ERF_Constants.H:61
constexpr amrex::Real PI
Definition: ERF_Constants.H:24
constexpr amrex::Real nzeror
Definition: ERF_Constants.H:76
constexpr amrex::Real rhos
Definition: ERF_Constants.H:48
constexpr amrex::Real rhor
Definition: ERF_Constants.H:47
constexpr amrex::Real a_rain
Definition: ERF_Constants.H:57
constexpr amrex::Real nzeros
Definition: ERF_Constants.H:77
constexpr amrex::Real a_snow
Definition: ERF_Constants.H:59
constexpr amrex::Real crain
Definition: ERF_Constants.H:100
Real rho_0
Definition: ERF_InitCustomPert_ABL.H:4
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMSurfaceAccumulation sam_surface_accumulation_from_component_fluxes(const SAMPrecipFluxComponents &component_fluxes, const amrex::Real dtn) noexcept
Definition: ERF_SAMUtils.H:1028
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int sam_precip_face_donor_k(const int face_k, const int k_lo, const int k_hi) noexcept
Definition: ERF_SAMUtils.H:809
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMPrecipFluxComponents sam_precip_flux_components_density_corrected(const SAMPrecipFluxComponents &precip_fluxes, const amrex::Real rho_0, const amrex::Real rho_avg) noexcept
Definition: ERF_SAMUtils.H:983
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMPrecipFluxComponents sam_precip_component_fluxes_from_face_state(const SAMPrecipComponentFaceState &face_state, const amrex::Real vrain, const amrex::Real vsnow, const amrex::Real vgrau) noexcept
Definition: ERF_SAMUtils.H:955
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE SAMPrecipComponentFaceState sam_precip_component_face_state(const amrex::Array4< const amrex::Real > &rho_array, const amrex::Array4< const amrex::Real > &tabs_array, const amrex::Array4< const amrex::Real > &qpr_array, const amrex::Array4< const amrex::Real > &qps_array, const amrex::Array4< const amrex::Real > &qpg_array, const int i, const int j, const int k, const int k_lo, const int k_hi) noexcept
Definition: ERF_SAMUtils.H:873
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real sam_limit_precip_component_flux(const amrex::Real raw_flux, const amrex::Real rho_donor, const amrex::Real q_donor, const amrex::Real detJ_donor, const amrex::Real coef) noexcept
Definition: ERF_SAMUtils.H:999
@ qp
Definition: ERF_Kessler.H:32
@ rain_accum
Definition: ERF_Kessler.H:34
@ qpg
Definition: ERF_Morrison.H:41
@ qps
Definition: ERF_Morrison.H:40
@ graup_accum
Definition: ERF_Morrison.H:52
@ qpr
Definition: ERF_Morrison.H:39
@ snow_accum
Definition: ERF_Morrison.H:51
Definition: ERF_SAMUtils.H:124
amrex::Real rho_avg
Definition: ERF_SAMUtils.H:125
Definition: ERF_SAMUtils.H:145
amrex::Real snow
Definition: ERF_SAMUtils.H:147
amrex::Real rain
Definition: ERF_SAMUtils.H:146
amrex::Real graupel
Definition: ERF_SAMUtils.H:148
Definition: ERF_SAMUtils.H:139
amrex::Real rain
Definition: ERF_SAMUtils.H:140
amrex::Real snow
Definition: ERF_SAMUtils.H:141
amrex::Real graupel
Definition: ERF_SAMUtils.H:142

Referenced by Advance().

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

◆ Qmoist_Ptr() [1/2]

const amrex::MultiFab* SAM::Qmoist_Ptr ( const int &  varIdx) const
inline
149  {
151  return mic_fab_vars[MicVarMap[varIdx]].get();
152  }
Here is the call graph for this function:

◆ Qmoist_Ptr() [2/2]

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

Reimplemented from NullMoist.

142  {
144  return mic_fab_vars[MicVarMap[varIdx]].get();
145  }
Here is the call graph for this function:

◆ Qmoist_Restart_Vars()

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

Reimplemented from NullMoist.

173  {
174  a_idx.clear();
175  a_names.clear();
176 
177  // NOTE: These are the indices to access into qmoist (not mic_fab_vars)
178  a_idx.push_back(0); a_names.push_back("RainAccum");
179  a_idx.push_back(1); a_names.push_back("SnowAccum");
180  a_idx.push_back(2); a_names.push_back("GraupAccum");
181  }

◆ Qmoist_Size()

int SAM::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

161 { return SAM::m_qmoist_size; }

◆ Qstate_Moist_NumConc_Size()

int SAM::Qstate_Moist_NumConc_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

int n_qstate_moist_numconc_size
Definition: ERF_SAM.H:306

◆ Qstate_Moist_Size()

int SAM::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

164 { return SAM::n_qstate_moist_size; }
int n_qstate_moist_size
Definition: ERF_SAM.H:303

◆ Set_dzmin()

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

Reimplemented from NullMoist.

101  {
102  m_dzmin = dz_min;
103  }

◆ Update_Micro_Vars()

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

Reimplemented from NullMoist.

115  {
116  this->Copy_State_to_Micro(cons_in);
117  this->Compute_Coefficients();
118  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitSAM.cpp:94
void Compute_Coefficients()
Definition: ERF_InitSAM.cpp:150
Here is the call graph for this function:

◆ Update_State_Vars()

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

Reimplemented from NullMoist.

123  {
124  this->Copy_Micro_to_State(cons_in);
125  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateSAM.cpp:15
Here is the call graph for this function:

Member Data Documentation

◆ accrgc

amrex::TableData<amrex::Real, 1> SAM::accrgc
private

◆ accrgi

amrex::TableData<amrex::Real, 1> SAM::accrgi
private

◆ accrrc

amrex::TableData<amrex::Real, 1> SAM::accrrc
private

◆ accrsc

amrex::TableData<amrex::Real, 1> SAM::accrsc
private

◆ accrsi

amrex::TableData<amrex::Real, 1> SAM::accrsi
private

◆ CFL_MAX

constexpr amrex::Real SAM::CFL_MAX = myhalf
staticconstexprprivate

◆ coefice

amrex::TableData<amrex::Real, 1> SAM::coefice
private

◆ dt

amrex::Real SAM::dt
private

Referenced by Advance().

◆ evapg1

amrex::TableData<amrex::Real, 1> SAM::evapg1
private

◆ evapg2

amrex::TableData<amrex::Real, 1> SAM::evapg2
private

◆ evapr1

amrex::TableData<amrex::Real, 1> SAM::evapr1
private

◆ evapr2

amrex::TableData<amrex::Real, 1> SAM::evapr2
private

◆ evaps1

amrex::TableData<amrex::Real, 1> SAM::evaps1
private

◆ evaps2

amrex::TableData<amrex::Real, 1> SAM::evaps2
private

◆ m_axis

int SAM::m_axis
private

Referenced by Define().

◆ m_detJ_cc

amrex::MultiFab* SAM::m_detJ_cc
private

◆ m_do_cond

bool SAM::m_do_cond
private

Referenced by Define().

◆ m_dzmin

amrex::Real SAM::m_dzmin
private

Referenced by Set_dzmin().

◆ m_fac_cond

amrex::Real SAM::m_fac_cond
private

Referenced by Define().

◆ m_fac_fus

amrex::Real SAM::m_fac_fus
private

Referenced by Define().

◆ m_fac_sub

amrex::Real SAM::m_fac_sub
private

Referenced by Define().

◆ m_geom

amrex::Geometry SAM::m_geom
private

◆ m_qmoist_size

int SAM::m_qmoist_size = 3
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_rdOcp

amrex::Real SAM::m_rdOcp
private

Referenced by Define().

◆ m_z_phys_nd

amrex::MultiFab* SAM::m_z_phys_nd
private

◆ mic_fab_vars

amrex::Array<FabPtr, MicVar::NumVars> SAM::mic_fab_vars
private

Referenced by Qmoist_Ptr().

◆ MicVarMap

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

Referenced by Qmoist_Ptr().

◆ n_qstate_moist_numconc_size

int SAM::n_qstate_moist_numconc_size = 0
private

◆ n_qstate_moist_size

int SAM::n_qstate_moist_size = 6
private

Referenced by Qstate_Moist_Size().

◆ nlev

int SAM::nlev
private

◆ pres1d

amrex::TableData<amrex::Real, 1> SAM::pres1d
private

◆ qn1d

amrex::TableData<amrex::Real, 1> SAM::qn1d
private

◆ qt1d

amrex::TableData<amrex::Real, 1> SAM::qt1d
private

◆ qv1d

amrex::TableData<amrex::Real, 1> SAM::qv1d
private

◆ rho1d

amrex::TableData<amrex::Real, 1> SAM::rho1d
private

◆ tabs1d

amrex::TableData<amrex::Real, 1> SAM::tabs1d
private

◆ zhi

int SAM::zhi
private

◆ zlo

int SAM::zlo
private

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