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

Functions

AMREX_GPU_HOST void compute_saturation_pressure_null (MultiFab &, const MultiFab &)
 
AMREX_GPU_HOST void compute_saturation_pressure_H2O (MultiFab &a_mf_sat_pressure, const MultiFab &a_mf_temperature)
 
AMREX_GPU_HOST void compute_saturation_vapfrac_null (MultiFab &, const MultiFab &)
 
AMREX_GPU_HOST void compute_saturation_vapfrac_H2O (MultiFab &a_mf_sat_vapfrac, const MultiFab &a_mf_temperature, const MultiFab &a_mf_pressure)
 
AMREX_GPU_HOST void compute_saturation_pressure_null (amrex::MultiFab &, const amrex::MultiFab &)
 Null function for saturation pressure. More...
 
AMREX_GPU_HOST void compute_saturation_pressure_H2O (amrex::MultiFab &, const amrex::MultiFab &)
 Compute saturation pressure for moist air. More...
 
AMREX_GPU_HOST void compute_saturation_vapfrac_null (amrex::MultiFab &, const amrex::MultiFab &, const amrex::MultiFab &)
 Null function for saturation vapour fraction. More...
 
AMREX_GPU_HOST void compute_saturation_vapfrac_H2O (amrex::MultiFab &, const amrex::MultiFab &, const amrex::MultiFab &)
 Compute saturation vapour fraction for moist air. More...
 

Detailed Description

Functions for computing saturation pressure and vapour fraction

Function Documentation

◆ compute_saturation_pressure_H2O() [1/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_pressure_H2O ( amrex::MultiFab &  ,
const amrex::MultiFab &   
)

Compute saturation pressure for moist air.

◆ compute_saturation_pressure_H2O() [2/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_pressure_H2O ( MultiFab &  a_mf_sat_pressure,
const MultiFab &  a_mf_temperature 
)
14  {
15  const auto& gvec = a_mf_sat_pressure.nGrowVect();
16  for (MFIter mfi(a_mf_sat_pressure, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
17  Box bx = mfi.tilebox();
18  bx.grow(gvec);
19  const Array4<Real>& psat_arr = a_mf_sat_pressure.array(mfi);
20  const Array4<Real const>& temperature_arr = a_mf_temperature.array(mfi);
21 
22  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
23  { psat_arr(i,j,k,0) = erf_esatw(temperature_arr(i,j,k,0))*100; } );
24  // formula gives pressure in hPa; we will save it in Pa.
25  }
26  }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+0.5) *dx[0];const Real y=(j+0.5) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - 1.0)/(2.0 *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(1.0+deltaT, inv_gm1);const Real T=(1.0+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=0.25 *(1.0+std::cos(PI *std::min(r2d_xy, R)/R));})
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real erf_esatw(amrex::Real t, bool use_empirical=false)
Definition: ERF_MicrophysicsUtils.H:68

Referenced by MaterialProperties::setProperties_H2O().

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

◆ compute_saturation_pressure_null() [1/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_pressure_null ( amrex::MultiFab &  ,
const amrex::MultiFab &   
)

Null function for saturation pressure.

◆ compute_saturation_pressure_null() [2/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_pressure_null ( MultiFab &  ,
const MultiFab &   
)
9 { }

◆ compute_saturation_vapfrac_H2O() [1/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_vapfrac_H2O ( amrex::MultiFab &  ,
const amrex::MultiFab &  ,
const amrex::MultiFab &   
)

Compute saturation vapour fraction for moist air.

◆ compute_saturation_vapfrac_H2O() [2/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_vapfrac_H2O ( MultiFab &  a_mf_sat_vapfrac,
const MultiFab &  a_mf_temperature,
const MultiFab &  a_mf_pressure 
)
35  {
36  const auto& gvec = a_mf_sat_vapfrac.nGrowVect();
37  for (MFIter mfi(a_mf_sat_vapfrac, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
38  Box bx = mfi.tilebox();
39  bx.grow(gvec);
40  const Array4<Real>& qsat_arr = a_mf_sat_vapfrac.array(mfi);
41  const Array4<Real const>& temperature_arr = a_mf_temperature.array(mfi);
42  const Array4<Real const>& pressure_arr = a_mf_pressure.array(mfi);
43 
44  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
45  {
46  // pressure is in Pa; formula takes pressure in hPa
47  erf_qsatw( temperature_arr(i,j,k,0),
48  pressure_arr(i,j,k,0)/Real(100.0),
49  qsat_arr(i,j,k,0) );
50  } );
51  }
52  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:171
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by MaterialProperties::setProperties_H2O().

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

◆ compute_saturation_vapfrac_null() [1/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_vapfrac_null ( amrex::MultiFab &  ,
const amrex::MultiFab &  ,
const amrex::MultiFab &   
)

Null function for saturation vapour fraction.

◆ compute_saturation_vapfrac_null() [2/2]

AMREX_GPU_HOST void saturation_funcs::compute_saturation_vapfrac_null ( MultiFab &  ,
const MultiFab &   
)
29 { }