ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_HSEUtils.H
Go to the documentation of this file.
1 #ifndef ERF_HSEUTIL_H_
2 #define ERF_HSEUTIL_H_
3 
5 #include "ERF_Constants.H"
6 #include "ERF_EOS.H"
7 
8 /**
9  * Utility functions for calculating a hydrostatic equilibrium (HSE) base state
10 */
11 
12 namespace HSEutils
13 {
14  using namespace amrex;
15 
16  // to control Newton iterations in init_isentropic_hse*
17  const int MAX_ITER = 10;
18  const amrex::Real TOL = 1.e-8;
19 
20  /**
21  * Function to calculate the hydrostatic density and pressure
22  * from Newton-Raphson iterations
23  *
24  * @param[in] m_tol iteration tolerance
25  * @param[in] RdOCp Rd/Cp
26  * @param[out] dz change in vertical height
27  * @param[out] g magnitude of gravity
28  * @param[in] C sum of known terms in HSE balance
29  * @param[in] T theta at the current cell center
30  * @param[in] qt total moisture (non-precip and precip)
31  * @param[in] qv water vapor
32  * @param[in] P pressure at cell center
33  * @param[in] rd dry density at cell center
34  * @param[in] F starting residual of non-linear eq
35  */
36  AMREX_GPU_HOST_DEVICE
37  AMREX_FORCE_INLINE
38  void
39  Newton_Raphson_hse (const Real& m_tol,
40  const Real& RdOCp,
41  const Real& dz,
42  const Real& g,
43  const Real& C,
44  const Real& T,
45  const Real& qt,
46  const Real& qv,
47  Real& P,
48  Real& rd,
49  Real& F)
50  {
51  int iter=0;
52  int max_iter=30;
53  Real eps = 1.0e-6;
54  Real ieps = 1.0e6;
55  do {
56  // Compute change in pressure
57  Real hi = getRhogivenThetaPress(T, P + eps, RdOCp, qv);
58  Real lo = getRhogivenThetaPress(T, P - eps, RdOCp, qv);
59  Real dFdp = 1.0 + 0.25*ieps*(hi - lo)*g*dz;
60  P -= F/dFdp;
61  if (P < 1.0e3) amrex::Warning("P < 1000 [Pa]; Domain height may be too large...");
62  AMREX_ALWAYS_ASSERT(P > 0.0);
63 
64  // Diagnose density and residual
65  rd = getRhogivenThetaPress(T, P, RdOCp, qv);
66  AMREX_ALWAYS_ASSERT(rd > 0.0);
67  Real r_tot = rd * (1. + qt);
68  F = P + 0.5*r_tot*g*dz + C;
69  ++iter;
70  }
71  while (std::abs(F)>m_tol && iter<max_iter);
72  if (iter>=max_iter) {
73  AMREX_DEVICE_PRINTF("WARNING: HSE Newton iterations did not converge to tolerance!\n%s", "");
74  AMREX_DEVICE_PRINTF("HSE Newton tol: %e %e\n",F,m_tol);
75  }
76  }
77 
78 
79  /**
80  * Function to calculate the hydrostatic density and pressure
81  *
82  * @param[in] r_sfc surface density
83  * @param[in] theta surface potential temperature
84  * @param[out] r hydrostatically balanced density profile
85  * @param[out] p hydrostatically balanced pressure profile
86  * @param[in] dz vertical grid spacing (constant)
87  * @param[in] prob_lo_z surface height
88  * @param[in] khi z-index corresponding to the big end of the domain
89  */
90  AMREX_GPU_HOST_DEVICE
91  AMREX_FORCE_INLINE
92  void
94  const amrex::Real& theta,
95  amrex::Real* r,
96  amrex::Real* p,
97  const amrex::Real& dz,
98  const int klo, const int khi)
99  {
100  int kstart;
101 
102  // r_sfc / p_0 are the density / pressure at klo
103 
104  // Initial guess
105  Real half_dz = 0.5*dz;
106  r[klo] = r_sfc;
107  p[klo] = p_0 - half_dz * r[klo] * CONST_GRAV;
108 
109  if (klo == 0)
110  {
111  kstart = 1;
112 
113  // We do a Newton iteration to satisfy the EOS & HSE (with constant theta)
114  bool converged_hse = false;
115  Real p_hse;
116  Real p_eos;
117 
118  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
119  {
120  p_hse = p_0 - half_dz * r[klo] * CONST_GRAV;
121  p_eos = getPgivenRTh(r[klo]*theta);
122 
123  Real A = p_hse - p_eos;
124 
126 
127  Real drho = A / (dpdr + half_dz * CONST_GRAV);
128 
129  r[klo] = r[klo] + drho;
130  p[klo] = getPgivenRTh(r[klo]*theta);
131 
132  if (std::abs(drho) < TOL)
133  {
134  converged_hse = true;
135  break;
136  }
137  }
138 
139  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << klo << std::endl;
140  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
141  } else {
142  kstart = klo; // because we need to solve for r[klo] here; r[klo-1] is what was passed in r_sfc
143  r[klo-1] = r_sfc; // r_sfc is really the interpolated density in the ghost cell in this case
144  p[klo-1] = getPgivenRTh(r[klo-1]*theta);
145  }
146 
147  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
148  for (int k = kstart; k <= khi; k++)
149  {
150  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
151  // to discretely satisfy HSE -- here we assume spatial_order = 2 -- we can generalize this later if needed
152  bool converged_hse = false;
153 
154  r[k] = r[k-1];
155 
156  Real p_eos = getPgivenRTh(r[k]*theta);
157  Real p_hse;
158 
159  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
160  {
161  Real r_avg = 0.5 * (r[k-1]+r[k]);
162  p_hse = p[k-1] - dz * r_avg * CONST_GRAV;
163  p_eos = getPgivenRTh(r[k]*theta);
164 
165  Real A = p_hse - p_eos;
166 
167  Real dpdr = getdPdRgivenConstantTheta(r[k],theta);
168  // Gamma * p_0 * std::pow( (R_d * theta / p_0), Gamma) * std::pow(r[k], Gamma-1.0) ;
169 
170  Real drho = A / (dpdr + dz * CONST_GRAV);
171 
172  r[k] = r[k] + drho;
173  p[k] = getPgivenRTh(r[k]*theta);
174 
175  if (std::abs(drho) < TOL * r[k-1])
176  {
177  converged_hse = true;
178  // amrex::Print() << " converged " << std::endl;
179  break;
180  }
181  }
182 
183  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k << std::endl;
184  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
185  }
186  r[khi+1] = r[khi];
187  }
188 
189 
190  /**
191  * Function to calculate the hydrostatic density and pressure over terrain
192  *
193  * @param[in] i x-index
194  * @param[in] j y-index
195  * @param[in] r_sfc surface density
196  * @param[in] theta surface potential temperature
197  * @param[out] r hydrostatically balanced density profile
198  * @param[out] p hydrostatically balanced pressure profile
199  * @param[in] z_cc cell-center heights
200  * @param[in] khi z-index corresponding to the big end of the domain
201  */
202  AMREX_GPU_HOST_DEVICE
203  AMREX_FORCE_INLINE
204  void
206  int j,
207  const amrex::Real& r_sfc,
208  const amrex::Real& theta,
209  amrex::Real* r,
210  amrex::Real* p,
211  const amrex::Array4<amrex::Real const> z_cc,
212  const int& klo, const int& khi)
213 {
214  int kstart;
215 
216  if (klo == 0) {
217  //
218  // r_sfc / p_0 are the density / pressure at the surface
219  //
220  // Initial guess
221  int k0 = 0;
222 
223  // Where we start the lower iteration
224  kstart = 1;
225 
226  Real half_dz = z_cc(i,j,k0);
227  r[k0] = r_sfc;
228  p[k0] = p_0 - half_dz * r[k0] * CONST_GRAV;
229  {
230  // We do a Newton iteration to satisfy the EOS & HSE (with constant theta)
231  bool converged_hse = false;
232  Real p_hse;
233  Real p_eos;
234 
235  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
236  {
237  p_hse = p_0 - half_dz * r[k0] * CONST_GRAV;
238  p_eos = getPgivenRTh(r[k0]*theta);
239 
240  Real A = p_hse - p_eos;
241 
242  Real dpdr = getdPdRgivenConstantTheta(r[k0],theta);
243 
244  Real drho = A / (dpdr + half_dz * CONST_GRAV);
245 
246  r[k0] = r[k0] + drho;
247  p[k0] = getPgivenRTh(r[k0]*theta);
248 
249  if (std::abs(drho) < TOL)
250  {
251  converged_hse = true;
252  break;
253  }
254  }
255 
256  //if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k0 << std::endl;
257  //if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
258  }
259  } else {
260  kstart = klo; // because we need to solve for r[klo] here; r[klo-1] is what was passed in r_sfc
261  r[klo-1] = r_sfc; // r_sfc is really the interpolated density in the ghost cell in this case
262  p[klo-1] = getPgivenRTh(r[klo-1]*theta);
263  }
264 
265  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
266  for (int k = kstart; k <= khi; k++)
267  {
268  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
269  // to discretely satisfy HSE -- here we assume spatial_order = 2 -- we can generalize this later if needed
270  bool converged_hse = false;
271 
272  Real dz_loc = (z_cc(i,j,k) - z_cc(i,j,k-1));
273 
274  r[k] = r[k-1];
275 
276  Real p_eos = getPgivenRTh(r[k]*theta);
277  Real p_hse;
278 
279  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
280  {
281  p_hse = p[k-1] - dz_loc * 0.5 * (r[k-1]+r[k]) * CONST_GRAV;
282  p_eos = getPgivenRTh(r[k]*theta);
283 
284  Real A = p_hse - p_eos;
285 
286  Real dpdr = getdPdRgivenConstantTheta(r[k],theta);
287 
288  Real drho = A / (dpdr + dz_loc * CONST_GRAV);
289 
290  r[k] = r[k] + drho;
291  p[k] = getPgivenRTh(r[k]*theta);
292 
293  if (std::abs(drho) < TOL * r[k-1])
294  {
295  converged_hse = true;
296  //amrex::Print() << " converged " << std::endl;
297  break;
298  }
299  }
300 
301  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k << std::endl;
302  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
303  }
304 }
305 
306 AMREX_FORCE_INLINE
307 AMREX_GPU_HOST_DEVICE
308 Real compute_saturation_pressure (const Real T_b, const bool use_empirical)
309 {
310  Real p_s = erf_esatw(T_b,use_empirical);
311  if (!use_empirical) p_s *= 100.0;
312  return p_s;
313 }
314 
315 AMREX_FORCE_INLINE
316 AMREX_GPU_HOST_DEVICE
317 Real compute_relative_humidity (const Real p_b, const Real T_b, const bool use_empirical,
318  const int which_zone, const Real scaled_height)
319 {
320  if (which_zone > 0) {
321  if (which_zone == 1) { // z <= height
322  Real p_s = compute_saturation_pressure(T_b, use_empirical);
323  Real q_s = Rd_on_Rv*p_s/(p_b - p_s);
324  return 0.014/q_s;
325  } else if (which_zone == 2) { // z > height and z <= z_tr; scaled_height = z/z_tr
326  return 1.0 - 0.75*std::pow(scaled_height,1.25);
327  } else { // z > z_tr
328  return 0.25;
329  }
330  } else {
331  return 1.0;
332  }
333 }
334 
335 AMREX_FORCE_INLINE
336 AMREX_GPU_HOST_DEVICE
337 Real vapor_mixing_ratio (const Real p_b, const Real T_b, const Real RH, const bool use_empirical,
338  int which_zone)
339 {
340  Real p_s = compute_saturation_pressure(T_b,use_empirical);
341  Real p_v = compute_vapor_pressure(p_s, RH);
342  Real q_v = Rd_on_Rv*p_v/(p_b - p_v);
343 
344  if (which_zone == 1) { // z <= height
345  return 0.014;
346  } else {
347  return q_v;
348  }
349 }
350 
351 AMREX_FORCE_INLINE
352 AMREX_GPU_HOST_DEVICE
353 Real compute_F_for_temp(const Real T_b, const Real p_b, const Real q_t, const Real eq_pot_temp,
354  const bool use_empirical, const int which_zone, const Real scaled_height)
355 {
356  Real fac = Cp_d + Cp_l*q_t;
357  Real RH = compute_relative_humidity(p_b, T_b, use_empirical, which_zone, scaled_height);
358  Real q_v = vapor_mixing_ratio(p_b, T_b, RH, use_empirical, which_zone);
359  Real p_s = compute_saturation_pressure(T_b, use_empirical);
360  Real p_v = compute_vapor_pressure(p_s, RH);
361  return eq_pot_temp - T_b*std::pow((p_b - p_v)/p_0, -R_d/fac)*std::exp(L_v*q_v/(fac*T_b));
362 }
363 
364 AMREX_FORCE_INLINE
365 AMREX_GPU_HOST_DEVICE
366 Real compute_temperature (const Real p_b, const Real q_t, const Real eq_pot_temp, const bool use_empirical,
367  const int which_zone, const Real scaled_height)
368 {
369  Real T_b = 200.0, delta_T; // Initial guess
370  for (int iter=0; iter<20; iter++)
371  {
372  Real F = compute_F_for_temp(T_b , p_b, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
373  Real F_plus_dF = compute_F_for_temp(T_b+1e-10, p_b, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
374  Real F_prime = (F_plus_dF - F)/1e-10;
375  delta_T = -F/F_prime;
376  T_b = T_b + delta_T;
377  }
378 
379  if(std::fabs(delta_T) > 1e-8){
380  amrex::Abort("Newton Raphson for temperature could not converge");
381  }
382 
383  return T_b;
384 }
385 
386 AMREX_FORCE_INLINE
387 AMREX_GPU_HOST_DEVICE
389 {
390  Real T_dp, gamma, T;
391  T = T_b - 273.15;
392 
393  Real b = 18.678, c = 257.14, d = 234.5;
394  gamma = log(RH*exp((b - T/d)*T/(c + T)));
395 
396  T_dp = c*gamma/(b - gamma);
397 
398  return T_dp;
399 }
400 
401 AMREX_FORCE_INLINE
402 AMREX_GPU_HOST_DEVICE
403 Real compute_theta (const Real scaled_height, const Real theta_0,
404  const Real theta_tr, const Real z_tr, const Real T_tr)
405 {
406  if(scaled_height <= 1.0) {
407  return theta_0 + (theta_tr - theta_0)*std::pow(scaled_height,1.25);
408  } else {
409  return theta_tr*exp(CONST_GRAV/(Cp_d*T_tr)*z_tr*(scaled_height - 1.0));
410  }
411 }
412 
413 AMREX_FORCE_INLINE
414 AMREX_GPU_HOST_DEVICE
415 void compute_rho (const Real& pressure, Real& theta, Real& rho, Real& q_v, Real& T_dp, Real& T_b,
416  const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone,
417  const Real scaled_height, const bool T_from_theta,
418  const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
419 {
420 
421  if (T_from_theta) {
422  theta = compute_theta(scaled_height, theta_0, theta_tr, z_tr, T_tr);
423  T_b = getTgivenPandTh(pressure, theta, (R_d/Cp_d));
424  } else {
425  T_b = compute_temperature(pressure, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
426  theta = getThgivenTandP(T_b, pressure, (R_d/Cp_d));
427  }
428 
429  Real RH = compute_relative_humidity(pressure, T_b, use_empirical, which_zone, scaled_height);
430 
431  q_v = vapor_mixing_ratio(pressure, T_b, RH, use_empirical, which_zone);
432 
433  rho = getRhogivenTandPress(T_b, pressure, q_v);
434 
435  if (T_from_theta) {
436  rho *= (1.0 + q_v);
437  } else {
438  rho *= (1.0 + q_t);
439  }
440 
441  T_dp = compute_dewpoint_temperature(T_b, RH);
442 }
443 
444 AMREX_FORCE_INLINE
445 AMREX_GPU_HOST_DEVICE
446 Real compute_F (const Real& p_k, const Real& p_k_minus_1, Real &theta_k, Real& rho_k, Real& q_v_k,
447  Real& T_dp, Real& T_b, const Real& dz, const Real& rho_k_minus_1, const Real q_t,
448  const Real eq_pot_temp, const bool use_empirical,
449  const int which_zone, const Real scaled_height, const bool T_from_theta,
450  const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
451 {
452  Real F;
453 
454  compute_rho(p_k, theta_k, rho_k, q_v_k, T_dp, T_b, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height,
455  T_from_theta, theta_0, theta_tr, z_tr, T_tr);
456 
457  if(rho_k_minus_1 == 0.0) // This loop is for the first point above the ground
458  {
459  F = p_k - p_k_minus_1 + rho_k*CONST_GRAV*dz/2.0;
460  }
461  else
462  {
463  F = p_k - p_k_minus_1 + 0.5 * (rho_k + rho_k_minus_1)*CONST_GRAV*dz;
464  }
465 
466  return F;
467 }
468 
469 AMREX_FORCE_INLINE
470 AMREX_GPU_HOST_DEVICE
471 Real compute_p_k (Real &p_k, const Real p_k_minus_1, Real &theta_k, Real& rho_k, Real& q_v_k,
472  Real& T_dp, Real& T_b, const Real dz, const Real rho_k_minus_1, const Real q_t,
473  const Real eq_pot_temp, const bool use_empirical,
474  const int which_zone, const Real scaled_height, const bool T_from_theta,
475  const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
476 {
477  Real delta_p_k;
478  Real eps = 1.e-10;
479 
480  for(int iter=0; iter<20; iter++)
481  {
482  Real F = compute_F(p_k , p_k_minus_1, theta_k, rho_k, q_v_k, T_dp, T_b, dz, rho_k_minus_1,
483  q_t, eq_pot_temp, use_empirical, which_zone, scaled_height, T_from_theta,
484  theta_0, theta_tr, z_tr, T_tr);
485  Real F_plus_dF = compute_F(p_k+eps, p_k_minus_1, theta_k, rho_k, q_v_k, T_dp, T_b, dz, rho_k_minus_1,
486  q_t, eq_pot_temp, use_empirical, which_zone, scaled_height, T_from_theta,
487  theta_0, theta_tr, z_tr, T_tr);
488  Real F_prime = (F_plus_dF - F)/eps;
489  delta_p_k = -F/F_prime;
490  p_k = p_k + delta_p_k;
491  }
492 
493 
494  if (std::fabs(delta_p_k) > 1e-8) {
495  amrex::Abort("Newton Raphson for pressure could not converge");
496  }
497 
498  return p_k;
499 }
500 
501 AMREX_FORCE_INLINE
502 AMREX_GPU_HOST_DEVICE
503 void
505  const Real& dz, const int& khi, const Real q_t,
506  const Real eq_pot_temp, const bool use_empirical,
507  const bool T_from_theta = false,
508  const Real z_tr_1 = -1., const Real z_tr_2 = -1.,
509  const Real theta_0 = 0., const Real theta_tr = 0., const Real T_tr = 0.)
510 {
511  Real T_b, T_dp;
512 
513  int which_zone = -1;
514  Real scaled_height = 0.;
515 
516  // Compute the quantities at z = 0.5*dz (first cell center)
517  p[0] = p_0;
518  if (z_tr_1 > 0.) {
519  Real z = 0.5 * dz;
520  if (z <= z_tr_1) {
521  which_zone = 1;
522  } else if (z <= z_tr_2) {
523  which_zone = 2;
524  } else {
525  which_zone = 3;
526  }
527  scaled_height = z/z_tr_2;
528  }
529 
530  compute_p_k(p[0], p_0, theta[0], r[0], q_v[0], T_dp, T_b, dz, 0.0, q_t, eq_pot_temp,
531  use_empirical, which_zone, scaled_height, T_from_theta,
532  theta_0, theta_tr, z_tr_2, T_tr);
533 
534  for (int k=1; k<=khi; k++)
535  {
536  p[k] = p[k-1];
537 
538  if (z_tr_1 > 0.) {
539  Real z = (k+0.5) * dz;
540  if (z <= z_tr_1) {
541  which_zone = 1;
542  } else if (z <= z_tr_2) {
543  which_zone = 2;
544  } else {
545  which_zone = 3;
546  }
547  scaled_height = z/z_tr_2;
548  }
549 
550  compute_p_k(p[k], p[k-1], theta[k], r[k], q_v[k], T_dp, T_b, dz, r[k-1], q_t, eq_pot_temp,
551  use_empirical, which_zone, scaled_height, T_from_theta,
552  theta_0, theta_tr, z_tr_2, T_tr);
553  }
554 
555  r[khi+1] = r[khi];
556 }
557 
558 } // namespace
559 
560 #endif
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real Rd_on_Rv
Definition: ERF_Constants.H:87
constexpr amrex::Real p_0
Definition: ERF_Constants.H:18
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
constexpr amrex::Real Cp_l
Definition: ERF_Constants.H:14
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
constexpr amrex::Real L_v
Definition: ERF_Constants.H:16
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 compute_vapor_pressure(const amrex::Real p_s, const amrex::Real RH)
Definition: ERF_EOS.H:181
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhogivenThetaPress(const amrex::Real th, const amrex::Real p, const amrex::Real rdOcp, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:96
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getdPdRgivenConstantTheta(const amrex::Real rho, const amrex::Real theta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:127
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhogivenTandPress(const amrex::Real T, const amrex::Real p, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:113
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenPandTh(const amrex::Real P, const amrex::Real th, const amrex::Real rdOcp)
Definition: ERF_EOS.H:32
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
amrex::Real gamma
Definition: ERF_InitCustomPert_IsentropicVortex.H:9
const int khi
Definition: ERF_InitCustomPert_ParticleTests.H:11
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real erf_esatw(amrex::Real t, bool use_empirical=false)
Definition: ERF_MicrophysicsUtils.H:68
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_HSEUtils.H:13
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void compute_rho(const Real &pressure, Real &theta, Real &rho, Real &q_v, Real &T_dp, Real &T_b, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone, const Real scaled_height, const bool T_from_theta, const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
Definition: ERF_HSEUtils.H:415
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_p_k(Real &p_k, const Real p_k_minus_1, Real &theta_k, Real &rho_k, Real &q_v_k, Real &T_dp, Real &T_b, const Real dz, const Real rho_k_minus_1, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone, const Real scaled_height, const bool T_from_theta, const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
Definition: ERF_HSEUtils.H:471
const int MAX_ITER
Definition: ERF_HSEUtils.H:17
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_dewpoint_temperature(const Real T_b, const Real RH)
Definition: ERF_HSEUtils.H:388
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real vapor_mixing_ratio(const Real p_b, const Real T_b, const Real RH, const bool use_empirical, int which_zone)
Definition: ERF_HSEUtils.H:337
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_saturation_pressure(const Real T_b, const bool use_empirical)
Definition: ERF_HSEUtils.H:308
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void Newton_Raphson_hse(const Real &m_tol, const Real &RdOCp, const Real &dz, const Real &g, const Real &C, const Real &T, const Real &qt, const Real &qv, Real &P, Real &rd, Real &F)
Definition: ERF_HSEUtils.H:39
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void init_isentropic_hse_no_terrain(Real *theta, Real *r, Real *p, Real *q_v, const Real &dz, const int &khi, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const bool T_from_theta=false, const Real z_tr_1=-1., const Real z_tr_2=-1., const Real theta_0=0., const Real theta_tr=0., const Real T_tr=0.)
Definition: ERF_HSEUtils.H:504
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_temperature(const Real p_b, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone, const Real scaled_height)
Definition: ERF_HSEUtils.H:366
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_relative_humidity(const Real p_b, const Real T_b, const bool use_empirical, const int which_zone, const Real scaled_height)
Definition: ERF_HSEUtils.H:317
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_F(const Real &p_k, const Real &p_k_minus_1, Real &theta_k, Real &rho_k, Real &q_v_k, Real &T_dp, Real &T_b, const Real &dz, const Real &rho_k_minus_1, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone, const Real scaled_height, const bool T_from_theta, const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
Definition: ERF_HSEUtils.H:446
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_F_for_temp(const Real T_b, const Real p_b, const Real q_t, const Real eq_pot_temp, const bool use_empirical, const int which_zone, const Real scaled_height)
Definition: ERF_HSEUtils.H:353
const amrex::Real TOL
Definition: ERF_HSEUtils.H:18
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void init_isentropic_hse_terrain(int i, int j, const amrex::Real &r_sfc, const amrex::Real &theta, amrex::Real *r, amrex::Real *p, const amrex::Array4< amrex::Real const > z_cc, const int &klo, const int &khi)
Definition: ERF_HSEUtils.H:205
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void init_isentropic_hse(const amrex::Real &r_sfc, const amrex::Real &theta, amrex::Real *r, amrex::Real *p, const amrex::Real &dz, const int klo, const int khi)
Definition: ERF_HSEUtils.H:93
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_theta(const Real scaled_height, const Real theta_0, const Real theta_tr, const Real z_tr, const Real T_tr)
Definition: ERF_HSEUtils.H:403
@ theta
Definition: ERF_MM5.H:20
@ rho
Definition: ERF_Kessler.H:22
@ qt
Definition: ERF_Kessler.H:27
@ qv
Definition: ERF_Kessler.H:28
@ T
Definition: ERF_IndexDefines.H:110
Definition: ERF_ConsoleIO.cpp:12
real(c_double), parameter g
Definition: ERF_module_model_constants.F90:19