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 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) 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_Size () override
 
void Qmoist_Restart_Vars (const SolverChoice &a_sc, std::vector< int > &a_idx, std::vector< std::string > &a_names) const override
 
- Public Member Functions inherited from NullMoist
 NullMoist ()
 
virtual ~NullMoist ()=default
 

Private Types

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

Private Attributes

int m_qmoist_size = 5
 
int m_qstate_size = 3
 
amrex::Vector< int > MicVarMap
 
amrex::Geometry m_geom
 
amrex::BoxArray m_gtoe
 
amrex::Real dt
 
int nlev
 
int zlo
 
int zhi
 
int m_axis
 
bool docloud
 
bool doprecip
 
amrex::Real m_fac_cond
 
amrex::Real m_fac_fus
 
amrex::Real m_fac_sub
 
amrex::Real m_gOcp
 
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 = 0.5
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ Kessler()

Kessler::Kessler ( )
inline
52 {}

◆ ~Kessler()

virtual Kessler::~Kessler ( )
virtualdefault

Member Function Documentation

◆ Advance()

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

Reimplemented from NullMoist.

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

◆ AdvanceKessler()

void Kessler::AdvanceKessler ( const SolverChoice solverChoice)

Compute Precipitation-related Microphysics quantities.

12 {
14  if (solverChoice.moisture_type == MoistureType::Kessler){
15  auto dz = m_geom.CellSize(2);
16  auto domain = m_geom.Domain();
17  int k_lo = domain.smallEnd(2);
18  int k_hi = domain.bigEnd(2);
19 
20  MultiFab fz;
21  auto ba = tabs->boxArray();
22  auto dm = tabs->DistributionMap();
23  fz.define(convert(ba, IntVect(0,0,1)), dm, 1, 0); // No ghost cells
24 
25  Real dtn = dt;
26  Real coef = dtn/dz;
27 
28  // Saturation and evaporation calculations go first
29  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
30  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
31  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
32  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
33  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
34  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
35  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
36  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
37  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
38 
39  const auto& tbx = mfi.tilebox();
40 
41  // Expose for GPU
42  Real d_fac_cond = m_fac_cond;
43 
44  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
45  {
46  qv_array(i,j,k) = std::max(0.0, qv_array(i,j,k));
47  qc_array(i,j,k) = std::max(0.0, qc_array(i,j,k));
48  qp_array(i,j,k) = std::max(0.0, qp_array(i,j,k));
49 
50  //------- Autoconversion/accretion
51  Real qcc, auto_r, accrr;
52  Real qsat, dtqsat;
53  Real dq_clwater_to_rain, dq_rain_to_vapor, dq_clwater_to_vapor, dq_vapor_to_clwater;
54 
55  Real pressure = pres_array(i,j,k);
56  erf_qsatw(tabs_array(i,j,k), pressure, qsat);
57  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat);
58 
59  if (qsat <= 0.0) {
60  amrex::Warning("qsat computed as non-positive; setting to 0.!");
61  qsat = 0.0;
62  }
63 
64  // If there is precipitating water (i.e. rain), and the cell is not saturated
65  // then the rain water can evaporate leading to extraction of latent heat, hence
66  // reducing temperature and creating negative buoyancy
67 
68  dq_clwater_to_rain = 0.0;
69  dq_rain_to_vapor = 0.0;
70  dq_vapor_to_clwater = 0.0;
71  dq_clwater_to_vapor = 0.0;
72 
73  //Real fac = qsat*4093.0*L_v/(Cp_d*std::pow(tabs_array(i,j,k)-36.0,2));
74  //Real fac = qsat*L_v*L_v/(Cp_d*R_v*tabs_array(i,j,k)*tabs_array(i,j,k));
75  Real fac = 1.0 + (L_v/Cp_d)*dtqsat;
76 
77  // If water vapor content exceeds saturation value, then vapor condenses to water and latent heat is released, increasing temperature
78  if (qv_array(i,j,k) > qsat) {
79  dq_vapor_to_clwater = std::min(qv_array(i,j,k), (qv_array(i,j,k)-qsat)/(1.0 + fac));
80  }
81 
82  // If water vapor is less than the saturated value, then the cloud water can evaporate,
83  // leading to evaporative cooling and reducing temperature
84  if (qv_array(i,j,k) < qsat && qc_array(i,j,k) > 0.0) {
85  dq_clwater_to_vapor = std::min(qc_array(i,j,k), (qsat - qv_array(i,j,k))/(1.0 + fac));
86  }
87 
88  if (qp_array(i,j,k) > 0.0 && qv_array(i,j,k) < qsat) {
89  Real C = 1.6 + 124.9*std::pow(0.001*rho_array(i,j,k)*qp_array(i,j,k),0.2046);
90  dq_rain_to_vapor = 1.0/(0.001*rho_array(i,j,k))*(1.0 - qv_array(i,j,k)/qsat)*C*std::pow(0.001*rho_array(i,j,k)*qp_array(i,j,k),0.525)/
91  (5.4e5 + 2.55e6/(pressure*qsat))*dtn;
92  // The negative sign is to make this variable (vapor formed from evaporation)
93  // a positive quantity (as qv/qs < 1)
94  dq_rain_to_vapor = std::min({qp_array(i,j,k), dq_rain_to_vapor});
95 
96  // Removing latent heat due to evaporation from rain water to water vapor, reduces the (potential) temperature
97  }
98 
99  // If there is cloud water present then do accretion and autoconversion to rain
100  if (qc_array(i,j,k) > 0.0) {
101  qcc = qc_array(i,j,k);
102 
103  auto_r = 0.0;
104  if (qcc > qcw0) {
105  auto_r = alphaelq;
106  }
107 
108  accrr = 0.0;
109  accrr = 2.2 * std::pow(qp_array(i,j,k) , 0.875);
110  dq_clwater_to_rain = dtn *(accrr*qcc + auto_r*(qcc - qcw0));
111 
112  // If the amount of change is more than the amount of qc present, then dq = qc
113  dq_clwater_to_rain = std::min(dq_clwater_to_rain, qc_array(i,j,k));
114  }
115 
116  qv_array(i,j,k) += -dq_vapor_to_clwater + dq_clwater_to_vapor + dq_rain_to_vapor;
117  qc_array(i,j,k) += dq_vapor_to_clwater - dq_clwater_to_vapor - dq_clwater_to_rain;
118  qp_array(i,j,k) += dq_clwater_to_rain - dq_rain_to_vapor;
119 
120  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
121  theta_array(i,j,k) += theta_over_T * d_fac_cond * (dq_vapor_to_clwater - dq_clwater_to_vapor - dq_rain_to_vapor);
122 
123  qv_array(i,j,k) = std::max(0.0, qv_array(i,j,k));
124  qc_array(i,j,k) = std::max(0.0, qc_array(i,j,k));
125  qp_array(i,j,k) = std::max(0.0, qp_array(i,j,k));
126 
127  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
128  });
129  }
130 
131  // Precompute terminal velocity for substepping
132  for ( MFIter mfi(fz, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
133  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
134  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
135 
136  auto fz_array = fz.array(mfi);
137  const Box& tbz = mfi.tilebox();
138 
139  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
140  {
141  Real rho_avg, qp_avg;
142 
143  if (k==k_lo) {
144  rho_avg = rho_array(i,j,k);
145  qp_avg = qp_array(i,j,k);
146  } else if (k==k_hi+1) {
147  rho_avg = rho_array(i,j,k-1);
148  qp_avg = qp_array(i,j,k-1);
149  } else {
150  rho_avg = 0.5*(rho_array(i,j,k-1) + rho_array(i,j,k)); // Convert to g/cm^3
151  qp_avg = 0.5*(qp_array(i,j,k-1) + qp_array(i,j,k));
152  }
153 
154  qp_avg = std::max(0.0, qp_avg);
155 
156  Real V_terminal = 36.34*std::pow(rho_avg*0.001*qp_avg, 0.1346)*std::pow(rho_avg/1.16, -0.5); // in m/s
157 
158  // NOTE: Fz is the sedimentation flux from the advective operator.
159  // In the terrain-following coordinate system, the z-deriv in
160  // the divergence uses the normal velocity (Omega). However,
161  // there are no u/v components to the sedimentation velocity.
162  // Therefore, we simply end up with a division by detJ when
163  // evaluating the source term: dJinv * (flux_hi - flux_lo) * dzinv.
164  fz_array(i,j,k) = rho_avg*V_terminal*qp_avg;
165  });
166  } // mfi
167 
168  // Compute number of substeps from maximum terminal velocity
169  Real wt_max;
170  int n_substep;
171  auto const& ma_fz_arr = fz.const_arrays();
172  GpuTuple<Real> max = ParReduce(TypeList<ReduceOpMax>{},
173  TypeList<Real>{},
174  fz, IntVect(0),
175  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
176  -> GpuTuple<Real>
177  {
178  return { ma_fz_arr[box_no](i,j,k) };
179  });
180  wt_max = get<0>(max) + std::numeric_limits<Real>::epsilon();
181  n_substep = int( std::ceil(wt_max * coef / CFL_MAX) );
182  AMREX_ALWAYS_ASSERT(n_substep >= 1);
183  coef /= Real(n_substep);
184  dtn /= Real(n_substep);
185 
186  // Substep the vertical advection
187  for (int nsub(0); nsub<n_substep; ++nsub) {
188  for ( MFIter mfi(*tabs, TilingIfNotGPU()); mfi.isValid(); ++mfi ){
189  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
190  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
191  auto rain_accum_array = mic_fab_vars[MicVar_Kess::rain_accum]->array(mfi);
192  auto fz_array = fz.array(mfi);
193 
194  const auto dJ_array = (m_detJ_cc) ? m_detJ_cc->const_array(mfi) : Array4<const Real>{};
195 
196  const Box& tbx = mfi.tilebox();
197  const Box& tbz = mfi.tilebox(IntVect(0,0,1),IntVect(0));
198 
199  // Update vertical flux every substep
200  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
201  {
202  Real rho_avg, qp_avg;
203  if (k==k_lo) {
204  rho_avg = rho_array(i,j,k);
205  qp_avg = qp_array(i,j,k);
206  } else if (k==k_hi+1) {
207  rho_avg = rho_array(i,j,k-1);
208  qp_avg = qp_array(i,j,k-1);
209  } else {
210  rho_avg = 0.5*(rho_array(i,j,k-1) + rho_array(i,j,k)); // Convert to g/cm^3
211  qp_avg = 0.5*(qp_array(i,j,k-1) + qp_array(i,j,k));
212  }
213 
214  qp_avg = std::max(0.0, qp_avg);
215 
216  Real V_terminal = 36.34*std::pow(rho_avg*0.001*qp_avg, 0.1346)*std::pow(rho_avg/1.16, -0.5); // in m/s
217 
218  // NOTE: Fz is the sedimentation flux from the advective operator.
219  // In the terrain-following coordinate system, the z-deriv in
220  // the divergence uses the normal velocity (Omega). However,
221  // there are no u/v components to the sedimentation velocity.
222  // Therefore, we simply end up with a division by detJ when
223  // evaluating the source term: dJinv * (flux_hi - flux_lo) * dzinv.
224  fz_array(i,j,k) = rho_avg*V_terminal*qp_avg;
225 
226  if(k==k_lo){
227  rain_accum_array(i,j,k) = rain_accum_array(i,j,k) + rho_avg*qp_avg*V_terminal*dtn/1000.0*1000.0; // Divide by rho_water and convert to mm
228  }
229  });
230 
231  // Update precip every substep
232  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
233  {
234  // Jacobian determinant
235  Real dJinv = (dJ_array) ? 1.0/dJ_array(i,j,k) : 1.0;
236 
237  if(std::fabs(fz_array(i,j,k+1)) < 1e-14) fz_array(i,j,k+1) = 0.0;
238  if(std::fabs(fz_array(i,j,k )) < 1e-14) fz_array(i,j,k ) = 0.0;
239  Real dq_sed = dJinv * (1.0/rho_array(i,j,k)) * (fz_array(i,j,k+1) - fz_array(i,j,k)) * coef;
240  if(std::fabs(dq_sed) < 1e-14) dq_sed = 0.0;
241 
242  qp_array(i,j,k) += dq_sed;
243  qp_array(i,j,k) = std::max(0.0, qp_array(i,j,k));
244  });
245  } // mfi
246  } // nsub
247  }
248 
249 
250  if (solverChoice.moisture_type == MoistureType::Kessler_NoRain){
251 
252  // get the temperature, dentisy, theta, qt and qc from input
253  for ( MFIter mfi(*tabs,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
254  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
255  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
256  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
257  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
258  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
259  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
260 
261  const auto& box3d = mfi.tilebox();
262 
263  // Expose for GPU
264  Real d_fac_cond = m_fac_cond;
265 
266  ParallelFor(box3d, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
267  {
268  qc_array(i,j,k) = std::max(0.0, qc_array(i,j,k));
269 
270  //------- Autoconversion/accretion
271  Real qsat, dtqsat;
272  Real dq_clwater_to_vapor, dq_vapor_to_clwater;
273 
274  Real pressure = pres_array(i,j,k);
275  erf_qsatw(tabs_array(i,j,k), pressure, qsat);
276  erf_dtqsatw(tabs_array(i,j,k), pressure, dtqsat);
277 
278  // If there is precipitating water (i.e. rain), and the cell is not saturated
279  // then the rain water can evaporate leading to extraction of latent heat, hence
280  // reducing temperature and creating negative buoyancy
281 
282  dq_vapor_to_clwater = 0.0;
283  dq_clwater_to_vapor = 0.0;
284 
285  //Real fac = qsat*4093.0*L_v/(Cp_d*std::pow(tabs_array(i,j,k)-36.0,2));
286  //Real fac = qsat*L_v*L_v/(Cp_d*R_v*tabs_array(i,j,k)*tabs_array(i,j,k));
287  Real fac = 1.0 + (L_v/Cp_d)*dtqsat;
288 
289  // If water vapor content exceeds saturation value, then vapor condenses to water and latent heat is released, increasing temperature
290  if (qv_array(i,j,k) > qsat){
291  dq_vapor_to_clwater = std::min(qv_array(i,j,k), (qv_array(i,j,k)-qsat)/(1.0 + fac));
292  }
293  // If water vapor is less than the saturated value, then the cloud water can evaporate, leading to evaporative cooling and
294  // reducing temperature
295  if (qv_array(i,j,k) < qsat && qc_array(i,j,k) > 0.0){
296  dq_clwater_to_vapor = std::min(qc_array(i,j,k), (qsat - qv_array(i,j,k))/(1.0 + fac));
297  }
298 
299  qv_array(i,j,k) += -dq_vapor_to_clwater + dq_clwater_to_vapor;
300  qc_array(i,j,k) += dq_vapor_to_clwater - dq_clwater_to_vapor;
301 
302  Real theta_over_T = theta_array(i,j,k)/tabs_array(i,j,k);
303 
304  theta_array(i,j,k) += theta_over_T * d_fac_cond * (dq_vapor_to_clwater - dq_clwater_to_vapor);
305 
306  qv_array(i,j,k) = std::max(0.0, qv_array(i,j,k));
307  qc_array(i,j,k) = std::max(0.0, qc_array(i,j,k));
308 
309  qt_array(i,j,k) = qv_array(i,j,k) + qc_array(i,j,k);
310  });
311  }
312  }
313 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real alphaelq
Definition: ERF_Constants.H:48
constexpr amrex::Real qcw0
Definition: ERF_Constants.H:46
constexpr amrex::Real L_v
Definition: ERF_Constants.H:16
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_dtqsatw(amrex::Real t, amrex::Real p, amrex::Real &dtqsatw)
Definition: ERF_MicrophysicsUtils.H:170
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:158
amrex::Geometry m_geom
Definition: ERF_Kessler.H:153
amrex::MultiFab * m_detJ_cc
Definition: ERF_Kessler.H:178
amrex::Real m_fac_cond
Definition: ERF_Kessler.H:171
amrex::Array< FabPtr, MicVar_Kess::NumVars > mic_fab_vars
Definition: ERF_Kessler.H:181
static constexpr amrex::Real CFL_MAX
Definition: ERF_Kessler.H:147
@ qp
Definition: ERF_Kessler.H:39
@ qcl
Definition: ERF_Kessler.H:37
@ tabs
Definition: ERF_Kessler.H:32
@ pres
Definition: ERF_Kessler.H:33
@ rho
Definition: ERF_Kessler.H:30
@ theta
Definition: ERF_Kessler.H:31
@ qt
Definition: ERF_Kessler.H:35
@ rain_accum
Definition: ERF_Kessler.H:41
@ qv
Definition: ERF_Kessler.H:36
MoistureType moisture_type
Definition: ERF_DataStruct.H:664

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
@ cons
Definition: ERF_IndexDefines.H:129

Referenced by Update_State_Vars().

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.

65 {
66  // Get the temperature, density, theta, qt and qp from input
67  for ( MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
68  const auto& box3d = mfi.tilebox();
69 
70  auto states_array = cons_in.array(mfi);
71 
72  auto qt_array = mic_fab_vars[MicVar_Kess::qt]->array(mfi);
73  auto qv_array = mic_fab_vars[MicVar_Kess::qv]->array(mfi);
74  auto qc_array = mic_fab_vars[MicVar_Kess::qcl]->array(mfi);
75 
76  auto qp_array = mic_fab_vars[MicVar_Kess::qp]->array(mfi);
77 
78  auto rho_array = mic_fab_vars[MicVar_Kess::rho]->array(mfi);
79  auto theta_array = mic_fab_vars[MicVar_Kess::theta]->array(mfi);
80  auto tabs_array = mic_fab_vars[MicVar_Kess::tabs]->array(mfi);
81  auto pres_array = mic_fab_vars[MicVar_Kess::pres]->array(mfi);
82 
83  // Get pressure, theta, temperature, density, and qt, qp
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)) * 0.01;
97  });
98  }
99 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:84
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#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.

63  {
64  docloud = sc.do_cloud;
65  doprecip = sc.do_precip;
66  m_fac_cond = lcond / sc.c_p;
67  m_fac_fus = lfus / sc.c_p;
68  m_fac_sub = lsub / sc.c_p;
69  m_gOcp = CONST_GRAV / sc.c_p;
70  m_axis = sc.ave_plane;
71  }
constexpr amrex::Real lfus
Definition: ERF_Constants.H:67
constexpr amrex::Real lsub
Definition: ERF_Constants.H:68
constexpr amrex::Real lcond
Definition: ERF_Constants.H:66
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
bool doprecip
Definition: ERF_Kessler.H:168
amrex::Real m_fac_sub
Definition: ERF_Kessler.H:173
bool docloud
Definition: ERF_Kessler.H:168
amrex::Real m_gOcp
Definition: ERF_Kessler.H:174
amrex::Real m_fac_fus
Definition: ERF_Kessler.H:172
int m_axis
Definition: ERF_Kessler.H:165
amrex::Real c_p
Definition: ERF_DataStruct.H:618
bool do_precip
Definition: ERF_DataStruct.H:678
bool do_cloud
Definition: ERF_DataStruct.H:677
int ave_plane
Definition: ERF_DataStruct.H:675

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

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

◆ Qmoist_Ptr()

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

Reimplemented from NullMoist.

116  {
117  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
118  return mic_fab_vars[MicVarMap[varIdx]].get();
119  }

◆ Qmoist_Restart_Vars()

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

Reimplemented from NullMoist.

131  {
132  a_idx.clear();
133  a_names.clear();
134  if (a_sc.moisture_type == MoistureType::Kessler) {
135  a_idx.push_back(4); a_names.push_back("RainAccum");
136  }
137  }

◆ Qmoist_Size()

int Kessler::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

122 { return Kessler::m_qmoist_size; }

◆ Qstate_Size()

int Kessler::Qstate_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

125 { return Kessler::m_qstate_size; }
int m_qstate_size
Definition: ERF_Kessler.H:144

◆ Update_Micro_Vars()

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

Reimplemented from NullMoist.

93  {
94  this->Copy_State_to_Micro(cons_in);
95  }
void Copy_State_to_Micro(const amrex::MultiFab &cons_in) override
Definition: ERF_InitKessler.cpp:64
Here is the call graph for this function:

◆ Update_State_Vars()

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

Reimplemented from NullMoist.

100  {
101  this->Copy_Micro_to_State(cons_in);
102  }
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 = 0.5
staticconstexprprivate

◆ docloud

bool Kessler::docloud
private

Referenced by Define().

◆ doprecip

bool Kessler::doprecip
private

Referenced by Define().

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

amrex::Real Kessler::m_fac_cond
private

Referenced by Define().

◆ m_fac_fus

amrex::Real Kessler::m_fac_fus
private

Referenced by Define().

◆ m_fac_sub

amrex::Real Kessler::m_fac_sub
private

Referenced by Define().

◆ m_geom

amrex::Geometry Kessler::m_geom
private

◆ m_gOcp

amrex::Real Kessler::m_gOcp
private

Referenced by Define().

◆ m_gtoe

amrex::BoxArray Kessler::m_gtoe
private

◆ m_qmoist_size

int Kessler::m_qmoist_size = 5
private

Referenced by Qmoist_Ptr(), and Qmoist_Size().

◆ m_qstate_size

int Kessler::m_qstate_size = 3
private

Referenced by Qstate_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().

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