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
 
void Set_RealWidth (const int real_width) override
 
void Qmoist_Restart_Vars (const SolverChoice &, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
 
- Public Member Functions inherited from NullMoist
 NullMoist ()
 
virtual ~NullMoist ()=default
 
virtual int Qstate_NonMoist_Size ()
 
virtual void GetPlotVarNames (amrex::Vector< std::string > &a_vec) const
 
virtual void GetPlotVar (const std::string &, amrex::MultiFab &) const
 
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
 

Private Types

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

Private Attributes

int m_qmoist_size = 1
 
int n_qstate_moist_size = 3
 
amrex::Vector< int > MicVarMap
 
amrex::Geometry m_geom
 
int m_real_width {0}
 
amrex::Real dt
 
amrex::Real m_dzmin
 
int nlev
 
int zlo
 
int zhi
 
int m_axis
 
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 = myhalf
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ Kessler()

Kessler::Kessler ( )
inline
44 {}

◆ ~Kessler()

virtual Kessler::~Kessler ( )
virtualdefault

Member Function Documentation

◆ Advance()

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

Reimplemented from NullMoist.

104  {
105  dt = dt_advance;
106 
107  this->AdvanceKessler(solverChoice);
108  }
amrex::Real dt
Definition: ERF_Kessler.H:158
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 i_lo = domain.smallEnd(0);
20  int i_hi = domain.bigEnd(0);
21  int j_lo = domain.smallEnd(1);
22  int j_hi = domain.bigEnd(1);
23  int k_lo = domain.smallEnd(2);
24  int k_hi = domain.bigEnd(2);
25  if (solverChoice.moisture_type == MoistureType::Kessler)
26  {
27  MultiFab fz;
28  auto ba = tabs->boxArray();
29  auto dm = tabs->DistributionMap();
30  fz.define(convert(ba, IntVect(0,0,1)), dm, 1, 0); // No ghost cells
31 
32  Real dtn = dt;
33  Real coef = dtn/m_dzmin;
34 
35  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
36  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
37  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
38  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
39  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
40  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
41  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
42  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
43  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
44 
45  auto tbx = mfi.tilebox();
46  if (tbx.smallEnd(0) == i_lo) { tbx.growLo(0,-m_real_width); }
47  if (tbx.bigEnd(0) == i_hi) { tbx.growHi(0,-m_real_width); }
48  if (tbx.smallEnd(1) == j_lo) { tbx.growLo(1,-m_real_width); }
49  if (tbx.bigEnd(1) == j_hi) { tbx.growHi(1,-m_real_width); }
50 
51  Real d_fac_cond = m_fac_cond;
52 
53  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
54  {
55  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
56  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
57  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
58 
59  Real qsat_local, dtqsat_local;
60  Real pressure = pres_array(i,j,k);
61  erf_qsatw(tabs_array(i,j,k), pressure, qsat_local);
62  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat_local);
63 
64  if (qsat_local <= Real(0)) {
65  amrex::Warning("qsat computed as non-positive; setting to Real(0)!");
66  qsat_local = Real(0);
67  }
68 
69  const KesslerSourceTerms source_terms = kessler_warm_rain_sources(
70  qv_array(i,j,k), qc_array(i,j,k), qp_array(i,j,k), rho_array(i,j,k),
71  pressure, qsat_local, dtqsat_local, dtn, do_cond);
72 
73  qv_array(i,j,k) += -source_terms.dq_vapor_to_cloud
74  + source_terms.dq_cloud_to_vapor
75  + source_terms.dq_rain_to_vapor;
76  qc_array(i,j,k) += source_terms.dq_vapor_to_cloud
77  - source_terms.dq_cloud_to_vapor
78  - source_terms.dq_cloud_to_rain;
79  qp_array(i,j,k) += source_terms.dq_cloud_to_rain
80  - source_terms.dq_rain_to_vapor;
81 
82  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
83  theta_array(i,j,k) += theta_over_T * d_fac_cond
84  * (source_terms.dq_vapor_to_cloud
85  - source_terms.dq_cloud_to_vapor
86  - source_terms.dq_rain_to_vapor);
87 
88  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
89  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
90  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
91 
92  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
93  });
94  }
95 
96  for ( MFIter mfi(fz, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
97  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
98  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
99 
100  auto fz_array = fz.array(mfi);
101  const Box& tbz = mfi.tilebox();
102 
103  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
104  {
105  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
106  const Real rho_k = (k == k_hi+1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
107  const Real qp_km1 = (k == k_lo) ? qp_array(i,j,k) : qp_array(i,j,k-1);
108  const Real qp_k = (k == k_hi+1) ? qp_array(i,j,k-1) : qp_array(i,j,k);
109  const KesslerFaceState face_state =
110  kessler_face_state(k, k_lo, k_hi, rho_km1, rho_k, qp_km1, qp_k);
111 
112  const Real terminal_velocity = kessler_terminal_velocity(face_state.rho, face_state.qp);
113  fz_array(i,j,k) = kessler_precip_flux(face_state.rho, terminal_velocity, face_state.qp);
114  });
115  }
116 
117  auto const& ma_fz_arr = fz.const_arrays();
118  GpuTuple<Real> max = ParReduce(TypeList<ReduceOpMax>{},
119  TypeList<Real>{},
120  fz, IntVect(0),
121  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
122  -> GpuTuple<Real>
123  {
124  return { ma_fz_arr[box_no](i,j,k) };
125  });
126  int n_substep = kessler_num_sedimentation_substeps(get<0>(max), dt, m_dzmin);
127  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(n_substep >= 1,
128  "Kessler: Number of precipitation substeps must be greater than 0!");
129  coef /= Real(n_substep);
130  dtn /= Real(n_substep);
131 
132  for (int nsub(0); nsub<n_substep; ++nsub) {
133  for ( MFIter mfi(*tabs, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
134  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
135  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
136  auto rain_accum_array = mic_fab_vars[MicVar_Kess::rain_accum]->array(mfi);
137  auto fz_array = fz.array(mfi);
138 
139  const auto dJ_array = (m_detJ_cc) ? m_detJ_cc->const_array(mfi) : Array4<const Real>{};
140 
141  const Box& tbx = mfi.tilebox();
142  const Box& tbz = mfi.tilebox(IntVect(0,0,1),IntVect(0));
143 
144  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
145  {
146  const Real rho_km1 = (k == k_lo) ? rho_array(i,j,k) : rho_array(i,j,k-1);
147  const Real rho_k = (k == k_hi+1) ? rho_array(i,j,k-1) : rho_array(i,j,k);
148  const Real qp_km1 = (k == k_lo) ? qp_array(i,j,k) : qp_array(i,j,k-1);
149  const Real qp_k = (k == k_hi+1) ? qp_array(i,j,k-1) : qp_array(i,j,k);
150  const KesslerFaceState face_state =
151  kessler_face_state(k, k_lo, k_hi, rho_km1, rho_k, qp_km1, qp_k);
152 
153  const Real terminal_velocity = kessler_terminal_velocity(face_state.rho, face_state.qp);
154  fz_array(i,j,k) = kessler_precip_flux(face_state.rho, terminal_velocity, face_state.qp);
155 
156  if(k==k_lo){
157  rain_accum_array(i,j,k) = rain_accum_array(i,j,k)
158  + face_state.rho * face_state.qp * terminal_velocity * dtn / Real(1000.0) * Real(1000.0);
159  }
160  });
161 
162  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
163  {
164  Real dJinv = (dJ_array) ? Real(1)/dJ_array(i,j,k) : Real(1);
165 
166  if (kessler_is_small_sedimentation_value(fz_array(i,j,k+1))) {
167  fz_array(i,j,k+1) = Real(0);
168  }
169  if (kessler_is_small_sedimentation_value(fz_array(i,j,k ))) {
170  fz_array(i,j,k ) = Real(0);
171  }
173  fz_array(i,j,k+1), fz_array(i,j,k), rho_array(i,j,k), dJinv, coef);
175  dq_sed = Real(0);
176  }
177 
178  qp_array(i,j,k) += dq_sed;
179  qp_array(i,j,k) = std::max(Real(0), qp_array(i,j,k));
180  });
181  }
182  }
183  }
184 
185  if (solverChoice.moisture_type == MoistureType::Kessler_NoRain) {
186  if (!do_cond) { return; }
187  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
188  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
189  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
190  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
191  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
192  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
193  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
194 
195  auto tbx = mfi.tilebox();
196  if (tbx.smallEnd(0) == i_lo) { tbx.growLo(0,-m_real_width); }
197  if (tbx.bigEnd(0) == i_hi) { tbx.growHi(0,-m_real_width); }
198  if (tbx.smallEnd(1) == j_lo) { tbx.growLo(1,-m_real_width); }
199  if (tbx.bigEnd(1) == j_hi) { tbx.growHi(1,-m_real_width); }
200 
201  Real d_fac_cond = m_fac_cond;
202 
203  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
204  {
205  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
206 
207  Real qsat, dtqsat;
208 
209  Real pressure = pres_array(i,j,k);
210  erf_qsatw(tabs_array(i,j,k), pressure, qsat);
211  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat);
212 
213  const KesslerSaturationAdjustment saturation_adjustment =
214  kessler_saturation_adjustment(qv_array(i,j,k), qc_array(i,j,k), qsat, dtqsat, do_cond);
215 
216  qv_array(i,j,k) += -saturation_adjustment.dq_vapor_to_cloud
217  + saturation_adjustment.dq_cloud_to_vapor;
218  qc_array(i,j,k) += saturation_adjustment.dq_vapor_to_cloud
219  - saturation_adjustment.dq_cloud_to_vapor;
220 
221  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
222  theta_array(i,j,k) += theta_over_T * d_fac_cond
223  * (saturation_adjustment.dq_vapor_to_cloud - saturation_adjustment.dq_cloud_to_vapor);
224 
225  qv_array(i,j,k) = std::max(Real(0), qv_array(i,j,k));
226  qc_array(i,j,k) = std::max(Real(0), qc_array(i,j,k));
227 
228  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
229  });
230  }
231  }
232 }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+myhalf) *dx[0];const Real y=(j+myhalf) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - one)/(two *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(one+deltaT, inv_gm1);const Real T=(one+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=fourth *(one+std::cos(PI *std::min(r2d_xy, R)/R));})
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) noexcept
Definition: ERF_KesslerUtils.H:67
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:50
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) noexcept
Definition: ERF_KesslerUtils.H:131
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE KesslerFaceState kessler_face_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 qp_km1, const amrex::Real qp_k) noexcept
Definition: ERF_KesslerUtils.H:92
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:120
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:32
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:42
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool kessler_is_small_sedimentation_value(const amrex::Real value) noexcept
Definition: ERF_KesslerUtils.H:60
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
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Geometry m_geom
Definition: ERF_Kessler.H:152
amrex::MultiFab * m_detJ_cc
Definition: ERF_Kessler.H:175
amrex::Real m_dzmin
Definition: ERF_Kessler.H:161
amrex::Real m_fac_cond
Definition: ERF_Kessler.H:170
bool m_do_cond
Definition: ERF_Kessler.H:171
amrex::Array< FabPtr, MicVar_Kess::NumVars > mic_fab_vars
Definition: ERF_Kessler.H:178
int m_real_width
Definition: ERF_Kessler.H:155
@ qp
Definition: ERF_Kessler.H:31
@ qcl
Definition: ERF_Kessler.H:29
@ tabs
Definition: ERF_Kessler.H:24
@ pres
Definition: ERF_Kessler.H:25
@ rho
Definition: ERF_Kessler.H:22
@ theta
Definition: ERF_Kessler.H:23
@ qt
Definition: ERF_Kessler.H:27
@ rain_accum
Definition: ERF_Kessler.H:33
@ qv
Definition: ERF_Kessler.H:28
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:1249

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
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
@ cons
Definition: ERF_IndexDefines.H:158

Referenced by Update_State_Vars().

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

◆ Copy_State_to_Micro()

void 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  // Get pressure, theta, temperature, density, and qt, qp
83  ParallelFor( box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
84  {
85  rho_array(i,j,k) = states_array(i,j,k,Rho_comp);
86  theta_array(i,j,k) = states_array(i,j,k,RhoTheta_comp)/states_array(i,j,k,Rho_comp);
87  qv_array(i,j,k) = states_array(i,j,k,RhoQ1_comp)/states_array(i,j,k,Rho_comp);
88  qc_array(i,j,k) = states_array(i,j,k,RhoQ2_comp)/states_array(i,j,k,Rho_comp);
89  qp_array(i,j,k) = states_array(i,j,k,RhoQ3_comp)/states_array(i,j,k,Rho_comp);
90  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
91 
92  tabs_array(i,j,k) = getTgivenRandRTh(states_array(i,j,k,Rho_comp),
93  states_array(i,j,k,RhoTheta_comp),
94  qv_array(i,j,k));
95  pres_array(i,j,k) = getPgivenRTh(states_array(i,j,k,RhoTheta_comp), qv_array(i,j,k)) * Real(0.01);
96  });
97  }
98 }
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.

55  {
56  m_fac_cond = lcond / sc.c_p;
57  m_axis = sc.ave_plane;
58  m_do_cond = (!sc.use_shoc);
59  }
constexpr amrex::Real lcond
Definition: ERF_Constants.H:77
int m_axis
Definition: ERF_Kessler.H:167
amrex::Real c_p
Definition: ERF_DataStruct.H:1168
bool use_shoc
Definition: ERF_DataStruct.H:1200
int ave_plane
Definition: ERF_DataStruct.H:1263

◆ 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:164
amrex::Vector< int > MicVarMap
Definition: ERF_Kessler.H:149
int zhi
Definition: ERF_Kessler.H:164
int nlev
Definition: ERF_Kessler.H:164
amrex::MultiFab * m_z_phys_nd
Definition: ERF_Kessler.H:174
int m_qmoist_size
Definition: ERF_Kessler.H:140
@ NumVars
Definition: ERF_Kessler.H:34
Here is the call graph for this function:

◆ Qmoist_Ptr()

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

Reimplemented from NullMoist.

112  {
114  return mic_fab_vars[MicVarMap[varIdx]].get();
115  }
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.

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

◆ Qmoist_Size()

int Kessler::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

118 { return Kessler::m_qmoist_size; }

◆ Qstate_Moist_Size()

int Kessler::Qstate_Moist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

int n_qstate_moist_size
Definition: ERF_Kessler.H:143

◆ Set_dzmin()

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

Reimplemented from NullMoist.

73  {
74  m_dzmin = dz_min;
75  }

◆ Set_RealWidth()

void Kessler::Set_RealWidth ( const int  real_width)
inlineoverridevirtual

Reimplemented from NullMoist.

124 { m_real_width = real_width; }

◆ Update_Micro_Vars()

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

96  {
97  this->Copy_Micro_to_State(cons_in);
98  }
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 = myhalf
staticconstexprprivate

◆ dt

amrex::Real Kessler::dt
private

Referenced by Advance().

◆ m_axis

int Kessler::m_axis
private

Referenced by Define().

◆ 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_real_width

int Kessler::m_real_width {0}
private

Referenced by Set_RealWidth().

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