ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
HSEutils Namespace Reference

Functions

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 &Th, const Real &qt, const Real &qv, Real &P, Real &rd, Real &F)
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void init_isentropic_hse_constant_dz (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)
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void init_isentropic_hse_stretched_dz (const amrex::Real &r_sfc, const amrex::Real &theta, amrex::Real *r, amrex::Real *p, const amrex::Real *stretched_dz, const int klo, const int khi)
 
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)
 
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_saturation_pressure (const Real T_b, const bool use_empirical)
 
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)
 
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)
 
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_F_for_temp_in_zone (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)
 
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)
 
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_dewpoint_temperature (const Real T_b, const Real RH)
 
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)
 
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)
 
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)
 
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_p_k_in_zone (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)
 
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=-one, const Real z_tr_2=-one, const Real theta_0=amrex::Real(0), const Real theta_tr=amrex::Real(0), const Real T_tr=amrex::Real(0))
 

Variables

const int MAX_ITER = 10
 
const amrex::Real TOL = Real(1.e-8)
 

Detailed Description

Utility functions for calculating a hydrostatic equilibrium (HSE) base state

Function Documentation

◆ compute_dewpoint_temperature()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_dewpoint_temperature ( const Real  T_b,
const Real  RH 
)
507 {
508  Real T_dp, gamma, T;
509  T = T_b - Real(273.15);
510 
511  Real b = Real(18.678), c = Real(257.14), d = Real(234.5);
512  gamma = std::log(RH*std::exp((b - T/d)*T/(c + T)));
513 
514  T_dp = c*gamma/(b - gamma);
515 
516  return T_dp;
517 }
Real T
Definition: ERF_InitCustomPert_Bubble.H:105
RH
Definition: ERF_InitCustomPert_Bubble.H:107
amrex::Real gamma
Definition: ERF_InitCustomPert_DataAssimilation_ISV.H:9
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by compute_rho().

Here is the caller graph for this function:

◆ compute_F()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::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 
)
569 {
570  Real F;
571 
572  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,
573  T_from_theta, theta_0, theta_tr, z_tr, T_tr);
574 
575  if(rho_k_minus_1 == amrex::Real(0)) // This loop is for the first point above the ground
576  {
577  F = p_k - p_k_minus_1 + rho_k*CONST_GRAV*dz/two;
578  }
579  else
580  {
581  F = p_k - p_k_minus_1 + myhalf * (rho_k + rho_k_minus_1)*CONST_GRAV*dz;
582  }
583 
584  return F;
585 }
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:40
Real eq_pot_temp
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:23
bool use_empirical
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:25
Real T_tr
Definition: ERF_InitCustomPert_SquallLine.H:43
Real q_t
Definition: ERF_InitCustomPert_SquallLine.H:24
Real theta_tr
Definition: ERF_InitCustomPert_SquallLine.H:47
Real theta_0
Definition: ERF_InitCustomPert_SquallLine.H:46
Real z_tr
Definition: ERF_InitCustomPert_SquallLine.H:36
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:533
@ dz
Definition: ERF_AdvanceWSM6.cpp:104

Referenced by compute_p_k_in_zone().

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

◆ compute_F_for_temp_in_zone()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_F_for_temp_in_zone ( 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 
)
467 {
468  Real fac = Cp_d + Cp_l*q_t;
469  Real RH = compute_relative_humidity(p_b, T_b, use_empirical, which_zone, scaled_height);
470  Real q_v = vapor_mixing_ratio(p_b, T_b, RH, use_empirical, which_zone);
472  Real p_v = compute_vapor_pressure(p_s, RH);
473  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));
474 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:31
constexpr amrex::Real p_0
Definition: ERF_Constants.H:37
constexpr amrex::Real Cp_l
Definition: ERF_Constants.H:33
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
constexpr amrex::Real L_v
Definition: ERF_Constants.H:35
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_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:449
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_saturation_pressure(const Real T_b, const bool use_empirical)
Definition: ERF_HSEUtils.H:421
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:429

Referenced by compute_temperature().

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

◆ compute_p_k_in_zone()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_p_k_in_zone ( 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 
)
594 {
595  Real delta_p_k;
596 
597 #ifdef AMREX_USE_FLOAT
598  Real eps = Real(1e-6);
599 #else
600  Real eps = Real(1e-10);
601 #endif
602 
603  for(int iter=0; iter<20; iter++)
604  {
605  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,
606  q_t, eq_pot_temp, use_empirical, which_zone, scaled_height, T_from_theta,
608  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,
609  q_t, eq_pot_temp, use_empirical, which_zone, scaled_height, T_from_theta,
611  Real F_prime = (F_plus_dF - F)/eps;
612  delta_p_k = -F/F_prime;
613  p_k = p_k + delta_p_k;
614  }
615 
616  if (std::fabs(delta_p_k) > TOL) {
617  amrex::Abort("Newton Raphson for pressure could not converge");
618  }
619 
620  return p_k;
621 }
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:564
const amrex::Real TOL
Definition: ERF_HSEUtils.H:21

Referenced by init_isentropic_hse_no_terrain().

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

◆ compute_relative_humidity()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_relative_humidity ( const Real  p_b,
const Real  T_b,
const bool  use_empirical,
const int  which_zone,
const Real  scaled_height 
)
431 {
432  if (which_zone > 0) {
433  if (which_zone == 1) { // z <= height
435  Real q_s = Rd_on_Rv*p_s/(p_b - p_s);
436  return Real(0.014)/q_s;
437  } else if (which_zone == 2) { // z > height and z <= z_tr; scaled_height = z/z_tr
438  return one - Real(0.75)*std::pow(scaled_height,Real(1.25));
439  } else { // z > z_tr
440  return fourth;
441  }
442  } else {
443  return one;
444  }
445 }
constexpr amrex::Real Rd_on_Rv
Definition: ERF_Constants.H:106
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12

Referenced by compute_F_for_temp_in_zone(), compute_rho(), and ParallelFor().

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

◆ compute_rho()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void HSEutils::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 
)
537 {
538 
539  if (T_from_theta) {
540  theta = compute_theta(scaled_height, theta_0, theta_tr, z_tr, T_tr);
541  T_b = getTgivenPandTh(pressure, theta, (R_d/Cp_d));
542  } else {
543  T_b = compute_temperature(pressure, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
544  theta = getThgivenTandP(T_b, pressure, (R_d/Cp_d));
545  }
546 
547  Real RH = compute_relative_humidity(pressure, T_b, use_empirical, which_zone, scaled_height);
548 
549  q_v = vapor_mixing_ratio(pressure, T_b, RH, use_empirical, which_zone);
550 
551  rho = getRhogivenTandPress(T_b, pressure, q_v);
552 
553  if (T_from_theta) {
554  rho *= (one + q_v);
555  } else {
556  rho *= (one + q_t);
557  }
558 
559  T_dp = compute_dewpoint_temperature(T_b, RH);
560 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhogivenTandPress(const amrex::Real T, const amrex::Real p, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:113
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 getTgivenPandTh(const amrex::Real P, const amrex::Real th, const amrex::Real rdOcp)
Definition: ERF_EOS.H:32
rho
Definition: ERF_InitCustomPert_Bubble.H:106
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_dewpoint_temperature(const Real T_b, const Real RH)
Definition: ERF_HSEUtils.H:506
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:478
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:521
@ theta
Definition: ERF_MM5.H:20

Referenced by compute_F().

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

◆ compute_saturation_pressure()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_saturation_pressure ( const Real  T_b,
const bool  use_empirical 
)
422 {
423  Real p_s = erf_esatw(T_b,use_empirical);
424  return p_s * Real(100.0);
425 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real erf_esatw(amrex::Real t, bool use_empirical=false)
Definition: ERF_MicrophysicsUtils.H:123

Referenced by compute_F_for_temp_in_zone(), compute_relative_humidity(), and vapor_mixing_ratio().

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

◆ compute_temperature()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::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 
)
480 {
481  Real T_b = Real(200.0), delta_T; // Initial guess
482 
483 #ifdef AMREX_USE_FLOAT
484  Real eps = Real(1.e-6) * T_b;
485 #else
486  Real eps = Real(1.e-10) * T_b;
487 #endif
488  for (int iter=0; iter<20; iter++)
489  {
490  Real F = compute_F_for_temp_in_zone(T_b , p_b, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
491  Real F_plus_dF = compute_F_for_temp_in_zone(T_b+eps, p_b, q_t, eq_pot_temp, use_empirical, which_zone, scaled_height);
492  Real F_prime = (F_plus_dF - F)/eps;
493  delta_T = -F/F_prime;
494  T_b = T_b + delta_T;
495  }
496 
497  if (std::fabs(delta_T) > TOL * T_b) {
498  amrex::Abort("Newton Raphson for temperature could not converge");
499  }
500 
501  return T_b;
502 }
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_F_for_temp_in_zone(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:465

Referenced by compute_rho().

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

◆ compute_theta()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::compute_theta ( const Real  scaled_height,
const Real  theta_0,
const Real  theta_tr,
const Real  z_tr,
const Real  T_tr 
)
523 {
524  if(scaled_height <= one) {
525  return theta_0 + (theta_tr - theta_0)*std::pow(scaled_height,Real(1.25));
526  } else {
527  return theta_tr*std::exp(CONST_GRAV/(Cp_d*T_tr)*z_tr*(scaled_height - one));
528  }
529 }

Referenced by compute_rho().

Here is the caller graph for this function:

◆ init_isentropic_hse_constant_dz()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void HSEutils::init_isentropic_hse_constant_dz ( 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 
)

Function to calculate the hydrostatic density and pressure with constant dz

Parameters
[in]r_sfcsurface density
[in]thetasurface potential temperature
[out]rhydrostatically balanced density profile
[out]phydrostatically balanced pressure profile
[in]dzvertical grid spacing (constant)
[in]kloz-index corresponding to the small end of the domain
[in]khiz-index corresponding to the big end of the domain
102  {
103  int kstart;
104 
105  // r_sfc / p_0 are the density / pressure at klo
106 
107  // Initial guess
108  Real myhalf_dz = myhalf*dz;
109  r[klo] = r_sfc;
110  p[klo] = p_0 - myhalf_dz * r[klo] * CONST_GRAV;
111 
112  if (klo == 0)
113  {
114  kstart = 1;
115 
116  // We do a Newton iteration to satisfy the EOS & HSE (with constant theta)
117  bool converged_hse = false;
118  Real p_hse;
119  Real p_eos;
120 
121  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
122  {
123  p_hse = p_0 - myhalf_dz * r[klo] * CONST_GRAV;
124  p_eos = getPgivenRTh(r[klo]*theta);
125 
126  Real A = p_hse - p_eos;
127 
128  Real dpdr = getdPdRgivenConstantTheta(r[klo],theta);
129 
130  Real drho = A / (dpdr + myhalf_dz * CONST_GRAV);
131 
132  r[klo] = r[klo] + drho;
133  p[klo] = getPgivenRTh(r[klo]*theta);
134 
135  if (std::abs(drho) < TOL)
136  {
137  converged_hse = true;
138  break;
139  }
140  }
141 
142  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << klo << std::endl;
143  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
144  } else {
145  kstart = klo; // because we need to solve for r[klo] here; r[klo-1] is what was passed in r_sfc
146  r[klo-1] = r_sfc; // r_sfc is really the interpolated density in the ghost cell in this case
147  p[klo-1] = getPgivenRTh(r[klo-1]*theta);
148  }
149 
150  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
151  for (int k = kstart; k <= khi; k++)
152  {
153  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
154  // to discretely satisfy HSE -- here we assume spatial_order = 2 -- we can generalize this later if needed
155  bool converged_hse = false;
156 
157  r[k] = r[k-1];
158 
159  Real p_eos = getPgivenRTh(r[k]*theta);
160  Real p_hse;
161 
162  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
163  {
164  Real r_avg = myhalf * (r[k-1]+r[k]);
165  p_hse = p[k-1] - dz * r_avg * CONST_GRAV;
166  p_eos = getPgivenRTh(r[k]*theta);
167 
168  Real A = p_hse - p_eos;
169 
170  Real dpdr = getdPdRgivenConstantTheta(r[k],theta);
171  // Gamma * p_0 * std::pow( (R_d * theta / p_0), Gamma) * std::pow(r[k], Gamma-one) ;
172 
173  Real drho = A / (dpdr + dz * CONST_GRAV);
174 
175  r[k] = r[k] + drho;
176  p[k] = getPgivenRTh(r[k]*theta);
177 
178  if (std::abs(drho) < TOL * r[k-1])
179  {
180  converged_hse = true;
181  // amrex::Print() << " converged " << std::endl;
182  break;
183  }
184  }
185 
186  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k << std::endl;
187  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
188  }
189  r[khi+1] = r[khi];
190  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getdPdRgivenConstantTheta(const amrex::Real rho, const amrex::Real theta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:127
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
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
const int MAX_ITER
Definition: ERF_HSEUtils.H:17
@ p
Definition: ERF_WSM6.H:176

Referenced by erf_init_dens_hse_dry().

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

◆ init_isentropic_hse_no_terrain()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void HSEutils::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 = -one,
const Real  z_tr_2 = -one,
const Real  theta_0 = amrex::Real(0),
const Real  theta_tr = amrex::Real(0),
const Real  T_tr = amrex::Real(0) 
)
632 {
633  Real T_b, T_dp;
634 
635  int which_zone = -1;
636  Real scaled_height = amrex::Real(0);
637 
638  // Compute the quantities at z = myhalf*dz (first cell center)
639  p[0] = p_0;
640  if (z_tr_1 > amrex::Real(0)) {
641  Real z = myhalf * dz;
642  if (z <= z_tr_1) {
643  which_zone = 1;
644  } else if (z <= z_tr_2) {
645  which_zone = 2;
646  } else {
647  which_zone = 3;
648  }
649  scaled_height = z/z_tr_2;
650  }
651 
652  compute_p_k_in_zone(p[0], p_0, theta[0], r[0], q_v[0], T_dp, T_b, dz, amrex::Real(0), q_t, eq_pot_temp,
653  use_empirical, which_zone, scaled_height, T_from_theta,
654  theta_0, theta_tr, z_tr_2, T_tr);
655 
656  for (int k=1; k<=khi; k++)
657  {
658  p[k] = p[k-1];
659 
660  if (z_tr_1 > amrex::Real(0)) {
661  Real z = (k+myhalf) * dz;
662  if (z <= z_tr_1) {
663  which_zone = 1;
664  } else if (z <= z_tr_2) {
665  which_zone = 2;
666  } else {
667  which_zone = 3;
668  }
669  scaled_height = z/z_tr_2;
670  }
671 
672  compute_p_k_in_zone(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,
673  use_empirical, which_zone, scaled_height, T_from_theta,
674  theta_0, theta_tr, z_tr_2, T_tr);
675  }
676 
677  r[khi+1] = r[khi];
678 }
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real compute_p_k_in_zone(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:589

Referenced by erf_init_dens_hse_moist(), and if().

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

◆ init_isentropic_hse_stretched_dz()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void HSEutils::init_isentropic_hse_stretched_dz ( const amrex::Real r_sfc,
const amrex::Real theta,
amrex::Real r,
amrex::Real p,
const amrex::Real stretched_dz,
const int  klo,
const int  khi 
)

Function to calculate the hydrostatic density and pressure with stretched dz

Parameters
[in]r_sfcsurface density
[in]thetasurface potential temperature
[out]rhydrostatically balanced density profile
[out]phydrostatically balanced pressure profile
[in]stretched_dzvertical grid spacing (stretched)
[in]kloz-index corresponding to the small end of the domain
[in]khiz-index corresponding to the big end of the domain
212  {
213  int kstart;
214 
215  // r_sfc / p_0 are the density / pressure at klo
216 
217  // Initial guess
218  r[klo] = r_sfc;
219  p[klo] = p_0 - myhalf*stretched_dz[klo] * r[klo] * CONST_GRAV;
220 
221  if (klo == 0)
222  {
223  kstart = 1;
224 
225  // We do a Newton iteration to satisfy the EOS & HSE (with constant theta)
226  bool converged_hse = false;
227  Real p_hse;
228  Real p_eos;
229 
230  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
231  {
232  p_hse = p_0 - myhalf*stretched_dz[klo] * r[klo] * CONST_GRAV;
233  p_eos = getPgivenRTh(r[klo]*theta);
234 
235  Real A = p_hse - p_eos;
236 
237  Real dpdr = getdPdRgivenConstantTheta(r[klo],theta);
238 
239  Real drho = A / (dpdr + myhalf*stretched_dz[klo] * CONST_GRAV);
240 
241  r[klo] = r[klo] + drho;
242  p[klo] = getPgivenRTh(r[klo]*theta);
243 
244  if (std::abs(drho) < TOL)
245  {
246  converged_hse = true;
247  break;
248  }
249  }
250 
251  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << klo << std::endl;
252  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
253  } else {
254  kstart = klo; // because we need to solve for r[klo] here; r[klo-1] is what was passed in r_sfc
255  r[klo-1] = r_sfc; // r_sfc is really the interpolated density in the ghost cell in this case
256  p[klo-1] = getPgivenRTh(r[klo-1]*theta);
257  }
258 
259  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
260  for (int k = kstart; k <= khi; k++)
261  {
262  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
263  // to discretely satisfy HSE -- here we assume spatial_order = 2 -- we can generalize this later if needed
264  bool converged_hse = false;
265 
266  r[k] = r[k-1];
267 
268  Real p_eos = getPgivenRTh(r[k]*theta);
269  Real p_hse;
270 
271  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
272  {
273  Real r_avg = myhalf * (r[k-1]+r[k]);
274  Real dz_avg = myhalf * (stretched_dz[k-1] + stretched_dz[k]);
275  p_hse = p[k-1] - dz_avg * r_avg * CONST_GRAV;
276  p_eos = getPgivenRTh(r[k]*theta);
277 
278  Real A = p_hse - p_eos;
279 
280  Real dpdr = getdPdRgivenConstantTheta(r[k],theta);
281  // Gamma * p_0 * std::pow( (R_d * theta / p_0), Gamma) * std::pow(r[k], Gamma-one) ;
282 
283  Real drho = A / (dpdr + dz_avg * CONST_GRAV);
284 
285  r[k] = r[k] + drho;
286  p[k] = getPgivenRTh(r[k]*theta);
287 
288  if (std::abs(drho) < TOL * r[k-1])
289  {
290  converged_hse = true;
291  // amrex::Print() << " converged " << std::endl;
292  break;
293  }
294  }
295 
296  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k << std::endl;
297  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
298  }
299  r[khi+1] = r[khi];
300  }

Referenced by erf_init_dens_hse_dry().

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

◆ init_isentropic_hse_terrain()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void HSEutils::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 
)

Function to calculate the hydrostatic density and pressure over terrain

Parameters
[in]ix-index
[in]jy-index
[in]r_sfcsurface density
[in]thetasurface potential temperature
[out]rhydrostatically balanced density profile
[out]phydrostatically balanced pressure profile
[in]z_cccell-center heights
[in]khiz-index corresponding to the big end of the domain
326  {
327  int kstart;
328 
329  if (klo == 0) {
330  //
331  // r_sfc / p_0 are the density / pressure at the surface
332  //
333  // Initial guess
334  int k0 = 0;
335 
336  // Where we start the lower iteration
337  kstart = 1;
338 
339  Real myhalf_dz = z_cc(i,j,k0);
340  r[k0] = r_sfc;
341  p[k0] = p_0 - myhalf_dz * r[k0] * CONST_GRAV;
342  {
343  // We do a Newton iteration to satisfy the EOS & HSE (with constant theta)
344  bool converged_hse = false;
345  Real p_hse;
346  Real p_eos;
347 
348  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
349  {
350  p_hse = p_0 - myhalf_dz * r[k0] * CONST_GRAV;
351  p_eos = getPgivenRTh(r[k0]*theta);
352 
353  Real A = p_hse - p_eos;
354 
355  Real dpdr = getdPdRgivenConstantTheta(r[k0],theta);
356 
357  Real drho = A / (dpdr + myhalf_dz * CONST_GRAV);
358 
359  r[k0] = r[k0] + drho;
360  p[k0] = getPgivenRTh(r[k0]*theta);
361 
362  if (std::abs(drho) < TOL)
363  {
364  converged_hse = true;
365  break;
366  }
367  }
368 
369  //if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k0 << std::endl;
370  //if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
371  }
372  } else {
373  kstart = klo; // because we need to solve for r[klo] here; r[klo-1] is what was passed in r_sfc
374  r[klo-1] = r_sfc; // r_sfc is really the interpolated density in the ghost cell in this case
375  p[klo-1] = getPgivenRTh(r[klo-1]*theta);
376  }
377 
378  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
379  for (int k = kstart; k <= khi; k++)
380  {
381  // To get values at k > 0 we do a Newton iteration to satisfy the EOS (with constant theta) and
382  // to discretely satisfy HSE -- here we assume spatial_order = 2 -- we can generalize this later if needed
383  bool converged_hse = false;
384 
385  Real dz_loc = (z_cc(i,j,k) - z_cc(i,j,k-1));
386 
387  r[k] = r[k-1];
388 
389  Real p_eos = getPgivenRTh(r[k]*theta);
390  Real p_hse;
391 
392  for (int iter = 0; iter < MAX_ITER && !converged_hse; iter++)
393  {
394  p_hse = p[k-1] - dz_loc * myhalf * (r[k-1]+r[k]) * CONST_GRAV;
395  p_eos = getPgivenRTh(r[k]*theta);
396 
397  Real A = p_hse - p_eos;
398 
399  Real dpdr = getdPdRgivenConstantTheta(r[k],theta);
400 
401  Real drho = A / (dpdr + dz_loc * CONST_GRAV);
402 
403  r[k] = r[k] + drho;
404  p[k] = getPgivenRTh(r[k]*theta);
405 
406  if (std::abs(drho) < TOL * r[k-1])
407  {
408  converged_hse = true;
409  //amrex::Print() << " converged " << std::endl;
410  break;
411  }
412  }
413 
414  // if (!converged_hse) amrex::Print() << "DOING ITERATIONS AT K = " << k << std::endl;
415  // if (!converged_hse) amrex::Error("Didn't converge the iterations in init");
416  }
417  }

Referenced by erf_init_dens_hse_dry().

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

◆ Newton_Raphson_hse()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void HSEutils::Newton_Raphson_hse ( const Real m_tol,
const Real RdOCp,
const Real dz,
const Real g,
const Real C,
const Real Th,
const Real qt,
const Real qv,
Real P,
Real rd,
Real F 
)

Function to calculate the hydrostatic density and pressure from Newton-Raphson iterations

Parameters
[in]m_toliteration tolerance
[in]RdOCpRd/Cp
[out]dzchange in vertical height
[out]gmagnitude of gravity
[in]Csum of known terms in HSE balance
[in]Ttheta at the current cell center
[in]qttotal moisture (non-precip and precip)
[in]qvwater vapor
[in]Ppressure at cell center
[in]rddry density at cell center
[in]Fstarting residual of non-linear eq
54  {
55  int iter=0;
56  int max_iter=20;
57  Real eps = Real(1.e-6);
58  Real ieps = Real(1.e6);
59  while (std::abs(F)>m_tol && iter<max_iter) {
60  // Compute change in pressure
61  Real hi = getRhogivenThetaPress(Th, P + eps, RdOCp, qv);
62  Real lo = getRhogivenThetaPress(Th, P - eps, RdOCp, qv);
63  Real dFdp = one + fourth*ieps*(hi - lo)*g*dz;
64  P -= F/dFdp;
65  if (P < Real(1.e3)) amrex::Warning("P < 1000 [Pa]; Domain height may be too large...");
67 
68  // Diagnose density and residual
69  rd = getRhogivenThetaPress(Th, P, RdOCp, qv);
71  Real r_tot = rd * (one + qt);
72  F = P + myhalf*r_tot*g*dz + C;
73  ++iter;
74  }
75  if (iter>=max_iter) {
76  AMREX_DEVICE_PRINTF("WARNING: HSE Newton iterations did not converge to tolerance!\n%s", "");
77  AMREX_DEVICE_PRINTF("HSE Newton tol: %e %e\n",F,m_tol);
78  }
79  }
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=amrex::Real(0))
Definition: ERF_EOS.H:96
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
@ P
Definition: ERF_IndexDefines.H:164
@ qt
Definition: ERF_Kessler.H:28
@ qv
Definition: ERF_Kessler.H:29
real(c_double), parameter g
Definition: ERF_module_model_constants.F90:19

Referenced by InputSoundingData::calc_rho_p().

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

◆ vapor_mixing_ratio()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE Real HSEutils::vapor_mixing_ratio ( const Real  p_b,
const Real  T_b,
const Real  RH,
const bool  use_empirical,
int  which_zone 
)
451 {
453  Real p_v = compute_vapor_pressure(p_s, RH);
454  Real q_v = Rd_on_Rv*p_v/(p_b - p_v);
455 
456  if (which_zone == 1) { // z <= height
457  return Real(0.014);
458  } else {
459  return q_v;
460  }
461 }

Referenced by compute_F_for_temp_in_zone(), compute_rho(), if(), and ParallelFor().

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

Variable Documentation

◆ MAX_ITER

◆ TOL