ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InitDensityHSE.H File Reference
#include "ERF_HSEUtils.H"
Include dependency graph for ERF_InitDensityHSE.H:

Go to the source code of this file.

Functions

void erf_init_dens_hse (amrex::MultiFab &rho_hse, std::unique_ptr< amrex::MultiFab > &z_phys_nd, std::unique_ptr< amrex::MultiFab > &z_phys_cc, amrex::Geometry const &geom) override
 
void erf_init_dens_hse_moist (amrex::MultiFab &rho_hse, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &geom) override
 
void erf_init_const_dens_hse (amrex::MultiFab &rho_hse) override
 
void erf_init_const_dens_and_th_hse (amrex::MultiFab &rho_hse, amrex::MultiFab &p_hse, amrex::MultiFab &pi_hse, amrex::MultiFab &th_hse, amrex::MultiFab &qv_hse, amrex::Real l_rdOcp) override
 

Function Documentation

◆ erf_init_const_dens_and_th_hse()

void erf_init_const_dens_and_th_hse ( amrex::MultiFab &  rho_hse,
amrex::MultiFab &  p_hse,
amrex::MultiFab &  pi_hse,
amrex::MultiFab &  th_hse,
amrex::MultiFab &  qv_hse,
amrex::Real  l_rdOcp 
)
override
198 {
199  amrex::Real rho_0 = base_parms.rho_0;
200  amrex::Real T_0 = base_parms.T_0;
201 
202  rho_hse.setVal(rho_0);
203  th_hse.setVal(T_0);
204 
205  amrex::Real rt0 = rho_0 * T_0;
206 
207  amrex::Real p0 = getPgivenRTh(rt0);
208  p_hse.setVal(p0);
209 
210  amrex::Real pi0 = getExnergivenRTh(rt0, l_rdOcp);
211 
212  pi_hse.setVal(pi0);
213 
214  // Default to zero background moisture
215  qv_hse.setVal(0.0);
216 }
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 getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:156
amrex::Real rho_0
Definition: ERF_InitCustomPert_IsentropicVortex.H:33
Real T_0
Definition: ERF_InitCustomPert_TaylorGreenVortex.H:4
amrex::Real Real
Definition: ERF_ShocInterface.H:19
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40
Here is the call graph for this function:

◆ erf_init_const_dens_hse()

void erf_init_const_dens_hse ( amrex::MultiFab &  rho_hse)
override
181 {
182  amrex::Real rho_0 = base_parms.rho_0;
183  for ( amrex::MFIter mfi(rho_hse, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
184  {
185  amrex::Array4<amrex::Real> rho_hse_arr = rho_hse.array(mfi);
186  const amrex::Box& gbx = mfi.growntilebox(1);
187  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
188  {
189  rho_hse_arr(i,j,k) = rho_0;
190  });
191  }
192 }
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
Here is the call graph for this function:

◆ erf_init_dens_hse()

void erf_init_dens_hse ( amrex::MultiFab &  rho_hse,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
std::unique_ptr< amrex::MultiFab > &  z_phys_cc,
amrex::Geometry const &  geom 
)
override

Initialize hydrostatically balanced density

Calls init_isentropic_hse_terrain() or init_isentropic_hse() for cases with and without terrain. Hydrostatic equilibrium (HSE) is satisfied discretely. Note that these routines presume that qv==0 when evaluating the EOS. Both density and pressure in HSE are calculated but only the density is used at this point.

18 {
19  const amrex::Real dz = geom.CellSize()[2];
20 
21  const amrex::Real T_sfc = 300.;
22  const amrex::Real rho_sfc = p_0 / (R_d*T_sfc);
23  const amrex::Real Thetabar = T_sfc;
24 
25  const int domlo_z = geom.Domain().smallEnd(2);
26  const int domhi_z = geom.Domain().bigEnd(2);
27 
28  if (domhi_z > 255) amrex::Abort("1D Arrays are hard-wired to only 256 high");
29 
30  bool use_terrain = (z_phys_nd != nullptr);
31 
32  if (use_terrain) {
33 
34  for ( amrex::MFIter mfi(rho_hse, TileNoZ()); mfi.isValid(); ++mfi )
35  {
36  amrex::Array4<amrex::Real > rho_arr = rho_hse.array(mfi);
37  amrex::Array4<amrex::Real const> z_cc_arr = z_phys_cc->const_array(mfi);
38 
39  // Create a flat box with same horizontal extent but only one cell in vertical
40  const amrex::Box& tbz = mfi.nodaltilebox(2);
41  amrex::Box b2d = tbz; // Copy constructor
42  b2d.grow(0,1); b2d.grow(1,1); // Grow by one in the lateral directions
43  b2d.setRange(2,0);
44 
45  const int ilo = tbz.smallEnd(0);
46  const int jlo = tbz.smallEnd(1);
47  const int klo = tbz.smallEnd(2);
48  const int khi = tbz.bigEnd(2)-1;
49 
50  amrex::Real rho_local_sfc;
51  if (klo == 0) {
52  rho_local_sfc = rho_sfc;
53  } else {
54  rho_local_sfc = rho_arr(ilo,jlo,klo-1);
55  }
56 
57  amrex::ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int)
58  {
59  amrex::Array1D<amrex::Real,0,255> r;
60  amrex::Array1D<amrex::Real,0,255> p;
61  HSEutils::init_isentropic_hse_terrain(i,j,rho_local_sfc,Thetabar,&(r(0)),&(p(0)),z_cc_arr,klo,khi);
62 
63  for (int k = klo; k <= khi; k++) {
64  rho_arr(i,j,k) = r(k);
65  }
66 
67  // Impose Neumann conditions at bottom and top of domain boundary
68  if (klo == domlo_z) {
69  rho_arr(i,j,domlo_z-1) = rho_arr(i,j,domlo_z);
70  }
71  if (khi == domhi_z) {
72  rho_arr(i,j,domhi_z+1) = rho_arr(i,j,domhi_z);
73  }
74  });
75  } // mfi
76 
77  } else {
78 
79  // Note that this integrates over the entire domain height,
80  // but only copies into each box as appropriate.
81  // There is no assumption that any one box spans the whole vertical height.
82  const int klo = geom.Domain().smallEnd(2);
83  const int khi = geom.Domain().bigEnd(2);
84 
85  // These are at cell centers (unstaggered)
86  amrex::Vector<amrex::Real> h_r(khi+2);
87  amrex::Vector<amrex::Real> h_p(khi+2);
88 
89  amrex::Gpu::DeviceVector<amrex::Real> d_r(khi+2);
90  amrex::Gpu::DeviceVector<amrex::Real> d_p(khi+2);
91 
92  HSEutils::init_isentropic_hse(rho_sfc,Thetabar,h_r.data(),h_p.data(),dz,klo,khi);
93 
94  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, h_r.begin(), h_r.end(), d_r.begin());
95  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, h_p.begin(), h_p.end(), d_p.begin());
96 
97  amrex::Real* r = d_r.data();
98 
99 #ifdef _OPENMP
100 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
101 #endif
102  for ( amrex::MFIter mfi(rho_hse, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
103  {
104  const amrex::Box& bx = mfi.growntilebox(1);
105  const amrex::Array4<amrex::Real> rho_hse_arr = rho_hse[mfi].array();
106  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
107  {
108  int kk = std::max(k,0);
109  rho_hse_arr(i,j,k) = r[kk];
110  });
111  } // mfi
112  } // no terrain
113 }
constexpr amrex::Real p_0
Definition: ERF_Constants.H:18
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
const bool use_terrain
Definition: ERF_InitCustomPertVels_Terrain3DHemisphere.H:26
const int khi
Definition: ERF_InitCustomPert_ParticleTests.H:11
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
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
Here is the call graph for this function:

◆ erf_init_dens_hse_moist()

void erf_init_dens_hse_moist ( amrex::MultiFab &  rho_hse,
std::unique_ptr< amrex::MultiFab > &  ,
amrex::Geometry const &  geom 
)
override
119 {
120  const amrex::Real dz = geom.CellSize()[2];
121  const int khi = geom.Domain().bigEnd()[2];
122 
123  // These are at cell centers (unstaggered)
124  amrex::Vector<amrex::Real> h_r(khi+2);
125  amrex::Vector<amrex::Real> h_p(khi+2);
126  amrex::Vector<amrex::Real> h_t(khi+2);
127  amrex::Vector<amrex::Real> h_q_v(khi+2);
128 
129  amrex::ParmParse pp("prob");
130  amrex::Real q_t = 0.02;
131  pp.query("qt_init",q_t);
132 
133  amrex::Real eq_pot_temp = 320.0;
134  pp.query("eq_pot_temp",eq_pot_temp);
135 
136  bool use_empirical = false;
137  pp.query("use_empirical_psat",use_empirical);
138 
139  amrex::Real z_tr_1 = -1.;
140  amrex::Real z_tr_2 = -1.;
141  pp.query("height", z_tr_1);
142  pp.query("z_tr", z_tr_2);
143 
144  // We only set them to zero to make sure they are initialized
145  amrex::Real T_tr = 0.0, theta_tr = 0.0, theta_0 = 0.0;
146 
147  bool T_from_theta = false;
148  pp.query("T_from_theta_in_moist_init", T_from_theta);
149 
150  if (T_from_theta) {
151  pp.get("theta_tr",theta_tr); pp.get("theta_0",theta_0);
152  pp.get("T_tr",T_tr);
153  }
154 
155  HSEutils::init_isentropic_hse_no_terrain( h_t.data(), h_r.data(), h_p.data(),
156  h_q_v.data(), dz, khi, q_t, eq_pot_temp, use_empirical,
157  T_from_theta, z_tr_1, z_tr_2,
158  theta_0, theta_tr, T_tr);
159 
160  amrex::Gpu::DeviceVector<amrex::Real> d_r(khi+2);
161  amrex::Gpu::copy(amrex::Gpu::hostToDevice, h_r.begin(), h_r.end(), d_r.begin());
162  amrex::Real* r = d_r.data();
163 
164 #ifdef _OPENMP
165 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
166 #endif
167  for ( amrex::MFIter mfi(rho_hse,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
168  {
169  const amrex::Box& bx = mfi.growntilebox(1);
170  const amrex::Array4<amrex::Real> rho_hse_arr = rho_hse[mfi].array();
171  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
172  {
173  int kk = std::max(k,0);
174  rho_hse_arr(i,j,k) = r[kk];
175  });
176  } // mfi
177 }
ParmParse pp("prob")
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
Here is the call graph for this function: