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) 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
 

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
 
amrex::Real m_fac_fus
 
amrex::Real m_fac_sub
 
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 = 0.5
 

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.

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

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

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.

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)) * 0.01;
96  });
97  }
98 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
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.

55  {
56  m_fac_cond = lcond / sc.c_p;
57  m_fac_fus = lfus / sc.c_p;
58  m_fac_sub = lsub / sc.c_p;
59  m_axis = sc.ave_plane;
60  m_do_cond = (!sc.use_shoc);
61  }
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
amrex::Real m_fac_sub
Definition: ERF_Kessler.H:173
amrex::Real m_fac_fus
Definition: ERF_Kessler.H:172
int m_axis
Definition: ERF_Kessler.H:168
amrex::Real c_p
Definition: ERF_DataStruct.H:1123
bool use_shoc
Definition: ERF_DataStruct.H:1155
int ave_plane
Definition: ERF_DataStruct.H:1186

◆ 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:165
amrex::Vector< int > MicVarMap
Definition: ERF_Kessler.H:150
int zhi
Definition: ERF_Kessler.H:165
int nlev
Definition: ERF_Kessler.H:165
amrex::MultiFab * m_z_phys_nd
Definition: ERF_Kessler.H:177
int m_qmoist_size
Definition: ERF_Kessler.H:141
@ 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.

113  {
114  AMREX_ALWAYS_ASSERT(varIdx < m_qmoist_size);
115  return mic_fab_vars[MicVarMap[varIdx]].get();
116  }

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

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

◆ Qmoist_Size()

int Kessler::Qmoist_Size ( )
inlineoverridevirtual

Reimplemented from NullMoist.

119 { 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:144

◆ Set_dzmin()

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

Reimplemented from NullMoist.

75  {
76  m_dzmin = dz_min;
77  }

◆ Set_RealWidth()

void Kessler::Set_RealWidth ( const int  real_width)
inlineoverridevirtual

Reimplemented from NullMoist.

125 { m_real_width = real_width; }

◆ 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)
inlineoverridevirtual

Reimplemented from NullMoist.

97  {
98  this->Copy_Micro_to_State(cons_in);
99  }
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

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