15 std::unique_ptr<amrex::MultiFab>& z_phys_nd,
16 std::unique_ptr<amrex::MultiFab>& z_phys_cc,
17 amrex::Geometry
const& geom)
override
25 const int domlo_z = geom.Domain().smallEnd(2);
26 const int domhi_z = geom.Domain().bigEnd(2);
28 if (domhi_z > 255) amrex::Abort(
"1D Arrays are hard-wired to only 256 high");
34 for ( amrex::MFIter mfi(rho_hse,
TileNoZ()); mfi.isValid(); ++mfi )
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);
40 const amrex::Box& tbz = mfi.nodaltilebox(2);
42 b2d.grow(0,1); b2d.grow(1,1);
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;
52 rho_local_sfc = rho_sfc;
54 rho_local_sfc = rho_arr(ilo,jlo,
klo-1);
59 amrex::Array1D<amrex::Real,0,255> r;
60 amrex::Array1D<amrex::Real,0,255> p;
63 for (
int k =
klo; k <=
khi; k++) {
64 rho_arr(i,j,k) = r(k);
69 rho_arr(i,j,domlo_z-1) = rho_arr(i,j,domlo_z);
72 rho_arr(i,j,domhi_z+1) = rho_arr(i,j,domhi_z);
82 const int klo = geom.Domain().smallEnd(2);
83 const int khi = geom.Domain().bigEnd(2);
86 amrex::Vector<amrex::Real> h_r(
khi+2);
87 amrex::Vector<amrex::Real> h_p(
khi+2);
89 amrex::Gpu::DeviceVector<amrex::Real> d_r(
khi+2);
90 amrex::Gpu::DeviceVector<amrex::Real> d_p(
khi+2);
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());
100 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
102 for ( amrex::MFIter mfi(rho_hse, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
104 const amrex::Box& bx = mfi.growntilebox(1);
105 const amrex::Array4<amrex::Real> rho_hse_arr = rho_hse[mfi].array();
108 int kk = std::max(k,0);
109 rho_hse_arr(i,j,k) = r[kk];
117 std::unique_ptr<amrex::MultiFab>& ,
118 amrex::Geometry
const& geom)
override
121 const int khi = geom.Domain().bigEnd()[2];
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);
129 amrex::ParmParse
pp(
"prob");
131 pp.query(
"qt_init",q_t);
134 pp.query(
"eq_pot_temp",eq_pot_temp);
136 bool use_empirical =
false;
137 pp.query(
"use_empirical_psat",use_empirical);
141 pp.query(
"height", z_tr_1);
142 pp.query(
"z_tr", z_tr_2);
145 amrex::Real T_tr = 0.0, theta_tr = 0.0, theta_0 = 0.0;
147 bool T_from_theta =
false;
148 pp.query(
"T_from_theta_in_moist_init", T_from_theta);
151 pp.get(
"theta_tr",theta_tr);
pp.get(
"theta_0",theta_0);
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);
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());
165 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
167 for ( amrex::MFIter mfi(rho_hse,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
169 const amrex::Box& bx = mfi.growntilebox(1);
170 const amrex::Array4<amrex::Real> rho_hse_arr = rho_hse[mfi].array();
173 int kk = std::max(k,0);
174 rho_hse_arr(i,j,k) = r[kk];
183 for ( amrex::MFIter mfi(rho_hse, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
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
189 rho_hse_arr(i,j,k) =
rho_0;
196 amrex::MultiFab& pi_hse, amrex::MultiFab& th_hse,
197 amrex::MultiFab& qv_hse,
amrex::Real l_rdOcp)
override
202 rho_hse.setVal(
rho_0);
constexpr amrex::Real p_0
Definition: ERF_Constants.H:18
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
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
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
const bool use_terrain
Definition: ERF_InitCustomPertVels_Terrain3DHemisphere.H:26
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);} } })
amrex::Real rho_0
Definition: ERF_InitCustomPert_IsentropicVortex.H:33
const int khi
Definition: ERF_InitCustomPert_ParticleTests.H:11
Real T_0
Definition: ERF_InitCustomPert_TaylorGreenVortex.H:4
void erf_init_const_dens_hse(amrex::MultiFab &rho_hse) override
Definition: ERF_InitDensityHSE.H:180
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
Definition: ERF_InitDensityHSE.H:14
void erf_init_dens_hse_moist(amrex::MultiFab &rho_hse, std::unique_ptr< amrex::MultiFab > &, amrex::Geometry const &geom) override
Definition: ERF_InitDensityHSE.H:116
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
Definition: ERF_InitDensityHSE.H:195
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
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_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
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40