ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InitDensityHSEDry.H File Reference

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
 

Function Documentation

◆ 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.

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