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

#include <ERF_Kessler.H>

Inheritance diagram for Kessler:
Collaboration diagram for Kessler:

Public Member Functions

 Kessler ()
 
virtual ~Kessler ()=default
 
void AdvanceKessler (const SolverChoice &solverChoice)
 
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 &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_RealWidth (const int)
 

Private Types

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

Private Attributes

int m_qmoist_size = 1
 
int n_qstate_moist_size = 3
 
int n_qstate_moist_numconc_size = 0
 
amrex::Vector< int > MicVarMap
 
amrex::Geometry m_geom
 
amrex::Real dt
 
amrex::Real m_dzmin
 
int nlev
 
int zlo
 
int zhi
 
amrex::Real m_fac_cond
 
bool m_do_cond
 
amrex::MultiFab * m_z_phys_nd
 
amrex::MultiFab * m_detJ_cc
 
amrex::Array< FabPtr, MicVar_Kess::NumVarsmic_fab_vars
 

Static Private Attributes

static constexpr amrex::Real CFL_MAX = kessler_sedimentation_cfl_max
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ Kessler()

Kessler::Kessler ( )
inline
45 {}

◆ ~Kessler()

virtual Kessler::~Kessler ( )
virtualdefault

Member Function Documentation

◆ Advance()

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

Reimplemented from NullMoist.

106  {
107  dt = dt_advance;
108 
109  this->AdvanceKessler(solverChoice);
110  }
amrex::Real dt
Definition: ERF_Kessler.H:160
void AdvanceKessler(const SolverChoice &solverChoice)
Definition: ERF_Kessler.cpp:14
Here is the call graph for this function:

◆ AdvanceKessler()

void Kessler::AdvanceKessler ( const SolverChoice solverChoice)

Compute Precipitation-related Microphysics quantities.

15 {
16  bool do_cond = m_do_cond;
18  auto domain = m_geom.Domain();
19  int k_lo = domain.smallEnd(2);
20  int k_hi = domain.bigEnd(2);
21  if (solverChoice.moisture_type == MoistureType::Kessler)
22  {
23  MultiFab fz;
24  auto ba = tabs->boxArray();
25  auto dm = tabs->DistributionMap();
26  fz.define(convert(ba, IntVect(0,0,1)), dm, 1, 0); // No ghost cells
27 
28  Real dtn = dt;
29  Real coef = dtn/m_dzmin;
30 
31  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
32  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
33  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
34  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
35  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
36  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
37  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
38  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
39  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
40 
41  auto tbx = mfi.tilebox();
42 
43  Real d_fac_cond = m_fac_cond;
44 
45  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
46  {
47  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
48  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
49  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
50 
51  Real qsat_local, dtqsat_local;
52  Real pressure = pres_array(i,j,k);
53  // Kessler stores pressure in mbar / hPa for the qsat helpers.
54  erf_qsatw(tabs_array(i,j,k), pressure, qsat_local);
55  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat_local);
56 
57  if (qsat_local <= Real(0)) {
58  amrex::Warning("qsat computed as non-positive; setting to Real(0)!");
59  qsat_local = Real(0);
60  }
61 
62  const KesslerSourceTerms source_terms = kessler_warm_rain_sources(
63  qv_array(i,j,k), qc_array(i,j,k), qp_array(i,j,k), rho_array(i,j,k),
64  pressure, qsat_local, dtqsat_local, dtn, do_cond, d_fac_cond);
65 
66  qv_array(i,j,k) += -source_terms.dq_vapor_to_cloud
67  + source_terms.dq_cloud_to_vapor
68  + source_terms.dq_rain_to_vapor;
69  qc_array(i,j,k) += source_terms.dq_vapor_to_cloud
70  - source_terms.dq_cloud_to_vapor
71  - source_terms.dq_cloud_to_rain;
72  qp_array(i,j,k) += source_terms.dq_cloud_to_rain
73  - source_terms.dq_rain_to_vapor;
74 
75  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
76  theta_array(i,j,k) += theta_over_T * d_fac_cond
77  * (source_terms.dq_vapor_to_cloud
78  - source_terms.dq_cloud_to_vapor
79  - source_terms.dq_rain_to_vapor);
80 
81  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
82  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
83  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
84 
85  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
86  });
87  }
88 
89  for ( MFIter mfi(fz, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
90  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
91  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
92 
93  auto fz_array = fz.array(mfi);
94  const Box& tbz = mfi.tilebox();
95 
96  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
97  {
98  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
99  const Real rho_k = (k == k_hi+1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
100  const Real qp_km1 = (k == k_lo) ? qp_array(i,j,k) : qp_array(i,j,k-1);
101  const Real qp_k = (k == k_hi+1) ? qp_array(i,j,k-1) : qp_array(i,j,k);
102  const KesslerFaceState face_state =
103  kessler_face_state(k, k_hi, rho_km1, rho_k, qp_km1, qp_k);
104 
105  const Real terminal_velocity = kessler_terminal_velocity(face_state.rho, face_state.qp);
106  fz_array(i,j,k) = kessler_precip_flux(face_state.rho, terminal_velocity, face_state.qp);
107  });
108  }
109 
110  auto const& ma_rho_arr = mic_fab_vars[MicVar_Kess::rho]->const_arrays();
111  auto const& ma_qp_arr = mic_fab_vars[MicVar_Kess::qp]->const_arrays();
112  // The sedimentation CFL uses fall speed, not precipitating mass flux.
113  // fz stores rho * Vt * qp for the flux divergence below, so the reduction
114  // recomputes Vt from the same face state used for the flux.
115  GpuTuple<Real> max_terminal_velocity = ParReduce(TypeList<ReduceOpMax>{},
116  TypeList<Real>{},
117  fz, IntVect(0),
118  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
119  -> GpuTuple<Real>
120  {
121  const auto& rho_arr = ma_rho_arr[box_no];
122  const auto& qp_arr = ma_qp_arr[box_no];
123  const Real rho_km1 = (k == k_lo) ? rho_arr(i,j,k) : rho_arr(i,j,k-1);
124  const Real rho_k = (k == k_hi+1) ? rho_arr(i,j,k-1) : rho_arr(i,j,k);
125  const Real qp_km1 = (k == k_lo) ? qp_arr(i,j,k) : qp_arr(i,j,k-1);
126  const Real qp_k = (k == k_hi+1) ? qp_arr(i,j,k-1) : qp_arr(i,j,k);
127  const KesslerFaceState face_state =
128  kessler_face_state(k, k_hi, rho_km1, rho_k, qp_km1, qp_k);
129  return { kessler_terminal_velocity(face_state.rho, face_state.qp) };
130  });
131  int n_substep = kessler_num_sedimentation_substeps(get<0>(max_terminal_velocity),
132  dt, m_dzmin);
133  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(n_substep >= 1,
134  "Kessler: Number of precipitation substeps must be greater than 0!");
135  coef /= Real(n_substep);
136  dtn /= Real(n_substep);
137 
138  for (int nsub(0); nsub<n_substep; ++nsub) {
139  for ( MFIter mfi(*tabs, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
140  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
141  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
142  auto rain_accum_array = mic_fab_vars[MicVar_Kess::rain_accum]->array(mfi);
143  auto fz_array = fz.array(mfi);
144 
145  const auto dJ_array = (m_detJ_cc) ? m_detJ_cc->const_array(mfi) : Array4<const Real>{};
146 
147  const Box& tbx = mfi.tilebox();
148  const Box& tbz = mfi.tilebox(IntVect(0,0,1),IntVect(0));
149 
150  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
151  {
152  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
153  const Real rho_k = (k == k_hi+1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
154  const Real qp_km1 = (k == k_lo) ? qp_array(i,j,k) : qp_array(i,j,k-1);
155  const Real qp_k = (k == k_hi+1) ? qp_array(i,j,k-1) : qp_array(i,j,k);
156  const int donor_k = kessler_face_donor_k(k, k_hi);
157  const KesslerFaceState face_state =
158  kessler_face_state(k, k_hi, rho_km1, rho_k, qp_km1, qp_k);
159 
160  const Real terminal_velocity = kessler_terminal_velocity(face_state.rho, face_state.qp);
161  const Real donor_rho = rho_array(i,j,donor_k);
162  const Real donor_qp = std::max(Real(0), qp_array(i,j,donor_k));
163  const Real donor_detJ = (dJ_array) ? dJ_array(i,j,donor_k) : Real(1);
164  // The face flux uses face rho and donor qp. The limiter uses donor rho,
165  // donor qp, and donor detJ because it caps how much rain can leave the donor
166  // cell in this substep.
167  // Cap outgoing flux by the donor cell's detJ-weighted available rain water:
168  // F * dt / dz <= rho_donor * qp_donor * detJ_donor.
169  // This keeps compressed cells from losing more rain than they contain.
170  const Real max_flux = donor_rho * donor_qp * donor_detJ / coef;
171  fz_array(i,j,k) = amrex::min(
172  kessler_precip_flux(face_state.rho, terminal_velocity, face_state.qp), max_flux);
173 
174  if(k==k_lo){
175  // Surface accumulation stores the bottom-face precipitation mass per area
176  // increment converted to liquid-water depth.
177  rain_accum_array(i,j,k) = rain_accum_array(i,j,k)
178  + kessler_rain_accumulation_increment(fz_array(i,j,k) * dtn);
179  }
180  });
181 
182  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
183  {
184  Real dJinv = (dJ_array) ? Real(1)/dJ_array(i,j,k) : Real(1);
185 
186  // Threshold local face-flux copies only. Neighboring cells share fz
187  // faces, so the cell update must not mutate fz while applying the
188  // small-value cutoff.
189  Real f_hi = fz_array(i,j,k+1);
190  Real f_lo = fz_array(i,j,k );
191 
193  f_hi = Real(0);
194  }
196  f_lo = Real(0);
197  }
199  f_hi, f_lo, rho_array(i,j,k), dJinv, coef);
201  dq_sed = Real(0);
202  }
203 
204  qp_array(i,j,k) += dq_sed;
205  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
206  });
207  }
208  }
209  }
210 
211  if (solverChoice.moisture_type == MoistureType::Kessler_NoRain) {
212  if (!do_cond) { return; }
213  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
214  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
215  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
216  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
217  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
218  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
219  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
220 
221  auto tbx = mfi.tilebox();
222 
223  Real d_fac_cond = m_fac_cond;
224 
225  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
226  {
227  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
228 
229  Real qsat, dtqsat;
230 
231  Real pressure = pres_array(i,j,k);
232  // Kessler stores pressure in mbar / hPa for the qsat helpers.
233  erf_qsatw(tabs_array(i,j,k), pressure, qsat);
234  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat);
235 
236  const KesslerSaturationAdjustment saturation_adjustment =
237  kessler_saturation_adjustment(qv_array(i,j,k), qc_array(i,j,k), qsat, dtqsat,
238  do_cond, d_fac_cond);
239 
240  qv_array(i,j,k) += -saturation_adjustment.dq_vapor_to_cloud
241  + saturation_adjustment.dq_cloud_to_vapor;
242  qc_array(i,j,k) += saturation_adjustment.dq_vapor_to_cloud
243  - saturation_adjustment.dq_cloud_to_vapor;
244 
245  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
246  theta_array(i,j,k) += theta_over_T * d_fac_cond
247  * (saturation_adjustment.dq_vapor_to_cloud - saturation_adjustment.dq_cloud_to_vapor);
248 
249  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
250  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
251 
252  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
253  });
254  }
255  }
256 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE KesslerSaturationAdjustment kessler_saturation_adjustment(const amrex::Real qv, const amrex::Real qc, const amrex::Real qsat, const amrex::Real dtqsat, const bool do_cond, const amrex::Real latent_over_cp) noexcept
Definition: ERF_KesslerUtils.H:88
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int kessler_num_sedimentation_substeps(const amrex::Real current_reduced_value, const amrex::Real dt, const amrex::Real dzmin) noexcept
Definition: ERF_KesslerUtils.H:61
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE KesslerFaceState kessler_face_state(const int k, const int k_hi, const amrex::Real rho_km1, const amrex::Real rho_k, const amrex::Real qp_km1, const amrex::Real qp_k) noexcept
Definition: ERF_KesslerUtils.H:113
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE KesslerSourceTerms kessler_warm_rain_sources(const amrex::Real qv, const amrex::Real qc, const amrex::Real qp, const amrex::Real rho, const amrex::Real pressure_current_units, const amrex::Real qsat, const amrex::Real dtqsat, const amrex::Real dt, const bool do_cond, const amrex::Real latent_over_cp) noexcept
Definition: ERF_KesslerUtils.H:155
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real kessler_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_KesslerUtils.H:144
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real kessler_terminal_velocity(const amrex::Real rho, const amrex::Real qp) noexcept
Definition: ERF_KesslerUtils.H:43
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real kessler_rain_accumulation_increment(const amrex::Real precip_mass_per_area) noexcept
Definition: ERF_KesslerUtils.H:71
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real kessler_precip_flux(const amrex::Real rho, const amrex::Real terminal_velocity, const amrex::Real qp) noexcept
Definition: ERF_KesslerUtils.H:53
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool kessler_is_small_sedimentation_value(const amrex::Real value) noexcept
Definition: ERF_KesslerUtils.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int kessler_face_donor_k(const int k, const int k_hi) noexcept
Definition: ERF_KesslerUtils.H:36
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_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:228
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::Real Real
Definition: ERF_ShocInterface.H:19
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
amrex::Geometry m_geom
Definition: ERF_Kessler.H:157
amrex::MultiFab * m_detJ_cc
Definition: ERF_Kessler.H:174
amrex::Real m_dzmin
Definition: ERF_Kessler.H:163
amrex::Real m_fac_cond
Definition: ERF_Kessler.H:169
bool m_do_cond
Definition: ERF_Kessler.H:170
amrex::Array< FabPtr, MicVar_Kess::NumVars > mic_fab_vars
Definition: ERF_Kessler.H:177
@ qp
Definition: ERF_Kessler.H:32
@ qcl
Definition: ERF_Kessler.H:30
@ tabs
Definition: ERF_Kessler.H:25
@ pres
Definition: ERF_Kessler.H:26
@ rho
Definition: ERF_Kessler.H:23
@ theta
Definition: ERF_Kessler.H:24
@ qt
Definition: ERF_Kessler.H:28
@ rain_accum
Definition: ERF_Kessler.H:34
@ qv
Definition: ERF_Kessler.H:29
Definition: ERF_KesslerUtils.H:26
amrex::Real rho
Definition: ERF_KesslerUtils.H:27
amrex::Real qp
Definition: ERF_KesslerUtils.H:28
Definition: ERF_KesslerUtils.H:21
amrex::Real dq_cloud_to_vapor
Definition: ERF_KesslerUtils.H:23
amrex::Real dq_vapor_to_cloud
Definition: ERF_KesslerUtils.H:22
Definition: ERF_KesslerUtils.H:14
amrex::Real dq_rain_to_vapor
Definition: ERF_KesslerUtils.H:18
amrex::Real dq_cloud_to_vapor
Definition: ERF_KesslerUtils.H:16
amrex::Real dq_cloud_to_rain
Definition: ERF_KesslerUtils.H:17
amrex::Real dq_vapor_to_cloud
Definition: ERF_KesslerUtils.H:15
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:

◆ Copy_Micro_to_State()

void Kessler::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.

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

Initializes the Microphysics module.

Parameters
[in]cons_inConserved variables input

Reimplemented from NullMoist.

64 {
65  // Get the temperature, density, theta, qt and qp from input
66  for ( MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
67  const auto& box3d = mfi.tilebox();
68 
69  auto states_array = cons_in.array(mfi);
70 
71  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
72  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
73  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
74 
75  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
76 
77  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
78  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
79  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
80  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
81 
82  // getPgivenRTh returns Pa. Kessler stores pressure in mbar / hPa for the
83  // qsat helper path, so convert here after forming temperature and density.
84  ParallelFor( box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
85  {
86  rho_array(i,j,k) = states_array(i,j,k,Rho_comp);
87  theta_array(i,j,k) = states_array(i,j,k,RhoTheta_comp)/states_array(i,j,k,Rho_comp);
88  qv_array(i,j,k) = states_array(i,j,k,RhoQ1_comp)/states_array(i,j,k,Rho_comp);
89  qc_array(i,j,k) = states_array(i,j,k,RhoQ2_comp)/states_array(i,j,k,Rho_comp);
90  qp_array(i,j,k) = states_array(i,j,k,RhoQ3_comp)/states_array(i,j,k,Rho_comp);
91  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
92 
93  tabs_array(i,j,k) = getTgivenRandRTh(states_array(i,j,k,Rho_comp),
94  states_array(i,j,k,RhoTheta_comp),
95  qv_array(i,j,k));
96  pres_array(i,j,k) = getPgivenRTh(states_array(i,j,k,RhoTheta_comp), qv_array(i,j,k)) * Real(0.01);
97  });
98  }
99 }
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 Kessler::Define ( SolverChoice sc)
inlineoverridevirtual

Reimplemented from NullMoist.

56  {
57  // Use one latent-over-cp factor for both the saturation solve and theta
58  // update. L_v comes from ERF_Constants.H; c_p comes from SolverChoice.
59  m_fac_cond = L_v / sc.c_p;
60  m_do_cond = (!sc.use_shoc);
61  }
constexpr amrex::Real L_v
Definition: ERF_Constants.H:35
amrex::Real c_p
Definition: ERF_DataStruct.H:1221
bool use_shoc
Definition: ERF_DataStruct.H:1253

◆ Init()

void Kessler::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.

27 {
28  dt = dt_advance;
29  m_geom = geom;
30 
31  m_z_phys_nd = z_phys_nd.get();
32  m_detJ_cc = detJ_cc.get();
33 
34  MicVarMap.resize(m_qmoist_size);
36 
37  // initialize microphysics variables
38  for (auto ivar = 0; ivar < MicVar_Kess::NumVars; ++ivar) {
39  mic_fab_vars[ivar] = std::make_shared<MultiFab>(cons_in.boxArray(), cons_in.DistributionMap(),
40  1, cons_in.nGrowVect());
41  mic_fab_vars[ivar]->setVal(0.);
42  }
43 
44  // Set class data members
45  for ( MFIter mfi(cons_in, TileNoZ()); mfi.isValid(); ++mfi) {
46  const auto& box3d = mfi.tilebox();
47 
48  const auto& lo = lbound(box3d);
49  const auto& hi = ubound(box3d);
50 
51  nlev = box3d.length(2);
52  zlo = lo.z;
53  zhi = hi.z;
54  }
55 }
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
int zlo
Definition: ERF_Kessler.H:166
amrex::Vector< int > MicVarMap
Definition: ERF_Kessler.H:154
int zhi
Definition: ERF_Kessler.H:166
int nlev
Definition: ERF_Kessler.H:166
amrex::MultiFab * m_z_phys_nd
Definition: ERF_Kessler.H:173
int m_qmoist_size
Definition: ERF_Kessler.H:142
@ NumVars
Definition: ERF_Kessler.H:35
Here is the call graph for this function:

◆ Qmoist_Ptr()

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

Reimplemented from NullMoist.

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

◆ Qmoist_Restart_Vars()

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

Reimplemented from NullMoist.

132  {
133  a_idx.clear();
134  a_names.clear();
135 
136  // NOTE: These are the indices to access into qmoist (not mic_fab_vars)
137  a_idx.push_back(0); a_names.push_back("RainAccum");
138  }

◆ Qmoist_Size()

int Kessler::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

120 { return Kessler::m_qmoist_size; }

◆ Qstate_Moist_NumConc_Size()

int Kessler::Qstate_Moist_NumConc_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

int n_qstate_moist_numconc_size
Definition: ERF_Kessler.H:148

◆ Qstate_Moist_Size()

int Kessler::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

int n_qstate_moist_size
Definition: ERF_Kessler.H:145

◆ Set_dzmin()

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

Reimplemented from NullMoist.

75  {
76  m_dzmin = dz_min;
77  }

◆ Update_Micro_Vars()

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

Reimplemented from NullMoist.

90  {
91  this->Copy_State_to_Micro(cons_in);
92  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitKessler.cpp:63
Here is the call graph for this function:

◆ Update_State_Vars()

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

Reimplemented from NullMoist.

98  {
99  this->Copy_Micro_to_State(cons_in);
100  }
void Copy_Micro_to_State(amrex::MultiFab &cons_in) override
Definition: ERF_UpdateKessler.cpp:14
Here is the call graph for this function:

Member Data Documentation

◆ CFL_MAX

constexpr amrex::Real Kessler::CFL_MAX = kessler_sedimentation_cfl_max
staticconstexprprivate

◆ dt

amrex::Real Kessler::dt
private

Referenced by Advance().

◆ m_detJ_cc

amrex::MultiFab* Kessler::m_detJ_cc
private

◆ m_do_cond

bool Kessler::m_do_cond
private

Referenced by Define().

◆ m_dzmin

amrex::Real Kessler::m_dzmin
private

Referenced by Set_dzmin().

◆ m_fac_cond

amrex::Real Kessler::m_fac_cond
private

Referenced by Define().

◆ m_geom

amrex::Geometry Kessler::m_geom
private

◆ m_qmoist_size

int Kessler::m_qmoist_size = 1
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_z_phys_nd

amrex::MultiFab* Kessler::m_z_phys_nd
private

◆ mic_fab_vars

amrex::Array<FabPtr, MicVar_Kess::NumVars> Kessler::mic_fab_vars
private

Referenced by Qmoist_Ptr().

◆ MicVarMap

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

Referenced by Qmoist_Ptr().

◆ n_qstate_moist_numconc_size

int Kessler::n_qstate_moist_numconc_size = 0
private

◆ n_qstate_moist_size

int Kessler::n_qstate_moist_size = 3
private

Referenced by Qstate_Moist_Size().

◆ nlev

int Kessler::nlev
private

◆ zhi

int Kessler::zhi
private

◆ zlo

int Kessler::zlo
private

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