ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
SHOCInterface Class Reference

#include <ERF_ShocInterface.H>

Collaboration diagram for SHOCInterface:

Classes

struct  Buffer
 
struct  SHOCPostprocess
 
struct  SHOCPreprocess
 

Public Member Functions

 SHOCInterface (const int &lev, SolverChoice &sc)
 
void set_grids (int &level, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons, amrex::MultiFab *xvel, amrex::MultiFab *yvel, amrex::MultiFab *zvel, amrex::Real *w_subsid, amrex::MultiFab *tau13, amrex::MultiFab *tau23, amrex::MultiFab *hfx3, amrex::MultiFab *qfx3, amrex::MultiFab *eddyDiffs, amrex::MultiFab *z_phys)
 
void initialize_impl ()
 
void run_impl (const Real dt)
 
void finalize_impl (const Real dt)
 
void alloc_buffers ()
 
void dealloc_buffers ()
 
void mf_to_kokkos_buffers ()
 
void kokkos_buffers_to_mf (const Real dt)
 
std::string name () const
 
void set_eddy_diffs ()
 
void set_diff_stresses ()
 
void add_fast_tend (amrex::Vector< amrex::MultiFab > &S_rhs)
 
void add_slow_tend (const amrex::MFIter &mfi, const amrex::Box &tbx, const amrex::Array4< amrex::Real > &cc_rhs_arr)
 

Protected Member Functions

void check_flux_state_consistency (const double dt)
 
void apply_turbulent_mountain_stress ()
 
void set_computed_group_impl ()
 
size_t requested_buffer_size_in_bytes () const
 
void init_buffers ()
 

Protected Attributes

amrex::Vector< int > m_col_offsets
 
Int m_num_cols = 0
 
Int m_num_layers = 0
 
Int m_npbl
 
Int m_nadv
 
Int m_num_tracers = 3
 
Int m_num_vel_comp = 2
 
Int hdtime
 
int m_lev
 
int m_step
 
amrex::Geometry m_geom
 
amrex::BoxArray m_ba
 
amrex::MultiFab * m_cons = nullptr
 
amrex::MultiFab * m_xvel = nullptr
 
amrex::MultiFab * m_yvel = nullptr
 
amrex::MultiFab * m_zvel = nullptr
 
amrex::Realm_w_subsid = nullptr
 
amrex::MultiFab * m_tau13 = nullptr
 
amrex::MultiFab * m_tau23 = nullptr
 
amrex::MultiFab * m_hfx3 = nullptr
 
amrex::MultiFab * m_qfx3 = nullptr
 
amrex::MultiFab * m_mu = nullptr
 
amrex::MultiFab * m_z_phys = nullptr
 
amrex::MultiFab c_tend
 
amrex::MultiFab u_tend
 
amrex::MultiFab v_tend
 
bool m_first_step = true
 
Buffer m_buffer
 
SHF::SHOCInput input
 
SHF::SHOCInputOutput input_output
 
SHF::SHOCOutput output
 
SHF::SHOCHistoryOutput history_output
 
SHF::SHOCRuntime runtime_options
 
SHOCPreprocess shoc_preprocess
 
SHOCPostprocess shoc_postprocess
 
ekat::WorkspaceManager< Spack, KT::Device > workspace_mgr
 
view_1d tot_buff_view
 
bool apply_tms = false
 
bool check_flux_state = false
 
bool extra_shoc_diags = false
 
bool column_conservation_check = false
 
view_2d omega
 
view_1d surf_sens_flux
 
sview_2d surf_mom_flux
 
view_1d surf_evap
 
view_2d T_mid
 
view_2d qv
 
view_1d surf_drag_coeff_tms
 
view_2d p_mid
 
view_2d p_int
 
view_2d pseudo_dens
 
view_1d phis
 
view_3d horiz_wind
 
view_2d sgs_buoy_flux
 
view_2d tk
 
view_2d cldfrac_liq
 
view_2d tke
 
view_2d qc
 
view_1d pblh
 
view_2d inv_qc_relvar
 
view_2d tkh
 
view_2d w_sec
 
view_2d cldfrac_liq_prev
 
view_1d ustar
 
view_1d obklen
 
view_2d brunt
 
view_2d shoc_mix
 
view_2d isotropy
 
view_2d shoc_cond
 
view_2d shoc_evap
 
view_2d wthl_sec
 
view_2d thl_sec
 
view_2d wqw_sec
 
view_2d qw_sec
 
view_2d uw_sec
 
view_2d vw_sec
 
view_2d w3
 
view_3d_strided qtracers
 
view_1d vapor_flux
 
view_1d water_flux
 
view_1d ice_flux
 
view_1d heat_flux
 

Private Types

using SHF = scream::shoc::Functions< Real, KokkosDefaultDevice >
 
using PF = scream::PhysicsFunctions< KokkosDefaultDevice >
 
using C = scream::physics::Constants< Real >
 
using KT = ekat::KokkosTypes< KokkosDefaultDevice >
 
using SC = scream::shoc::Constants< Real >
 
using Spack = typename SHF::Spack
 
using IntSmallPack = typename SHF::IntSmallPack
 
using Smask = typename SHF::Smask
 
using view_1d_int = typename KT::template view_1d< Int >
 
using view_1d = typename SHF::view_1d< Real >
 
using view_1d_const = typename SHF::view_1d< const Real >
 
using view_2d = typename SHF::view_2d< SHF::Spack >
 
using view_2d_const = typename SHF::view_2d< const Spack >
 
using sview_2d = typename KT::template view_2d< Real >
 
using sview_2d_const = typename KT::template view_2d< const Real >
 
using view_3d = typename SHF::view_3d< Spack >
 
using view_3d_const = typename SHF::view_3d< const Spack >
 
using view_3d_strided = typename SHF::view_3d_strided< Spack >
 
using WSM = ekat::WorkspaceManager< Spack, KT::Device >
 
template<typename ScalarT >
using uview_1d = ekat::Unmanaged< typename KT::template view_1d< ScalarT > >
 
template<typename ScalarT >
using uview_2d = ekat::Unmanaged< typename KT::template view_2d< ScalarT > >
 

Member Typedef Documentation

◆ C

using SHOCInterface::C = scream::physics::Constants<Real>
private

◆ IntSmallPack

using SHOCInterface::IntSmallPack = typename SHF::IntSmallPack
private

◆ KT

using SHOCInterface::KT = ekat::KokkosTypes<KokkosDefaultDevice>
private

◆ PF

using SHOCInterface::PF = scream::PhysicsFunctions<KokkosDefaultDevice>
private

◆ SC

using SHOCInterface::SC = scream::shoc::Constants<Real>
private

◆ SHF

using SHOCInterface::SHF = scream::shoc::Functions<Real, KokkosDefaultDevice>
private

◆ Smask

using SHOCInterface::Smask = typename SHF::Smask
private

◆ Spack

using SHOCInterface::Spack = typename SHF::Spack
private

◆ sview_2d

using SHOCInterface::sview_2d = typename KT::template view_2d<Real>
private

◆ sview_2d_const

using SHOCInterface::sview_2d_const = typename KT::template view_2d<const Real>
private

◆ uview_1d

template<typename ScalarT >
using SHOCInterface::uview_1d = ekat::Unmanaged<typename KT::template view_1d<ScalarT> >
private

◆ uview_2d

template<typename ScalarT >
using SHOCInterface::uview_2d = ekat::Unmanaged<typename KT::template view_2d<ScalarT> >
private

◆ view_1d

using SHOCInterface::view_1d = typename SHF::view_1d<Real>
private

◆ view_1d_const

using SHOCInterface::view_1d_const = typename SHF::view_1d<const Real>
private

◆ view_1d_int

using SHOCInterface::view_1d_int = typename KT::template view_1d<Int>
private

◆ view_2d

using SHOCInterface::view_2d = typename SHF::view_2d<SHF::Spack>
private

◆ view_2d_const

using SHOCInterface::view_2d_const = typename SHF::view_2d<const Spack>
private

◆ view_3d

using SHOCInterface::view_3d = typename SHF::view_3d<Spack>
private

◆ view_3d_const

using SHOCInterface::view_3d_const = typename SHF::view_3d<const Spack>
private

◆ view_3d_strided

using SHOCInterface::view_3d_strided = typename SHF::view_3d_strided<Spack>
private

◆ WSM

using SHOCInterface::WSM = ekat::WorkspaceManager<Spack, KT::Device>
private

Constructor & Destructor Documentation

◆ SHOCInterface()

SHOCInterface::SHOCInterface ( const int &  lev,
SolverChoice sc 
)
8 {
9  //
10  // Defaults set from E3SM/components/eamxx/cime_config/namelist_defaults_eamxx.xml
11  //
12  // Turn off SGS variability in SHOC, effectively reducing it to a Real(1.5) TKE closure?
13  bool def_shoc_1p5tke = false;
14  // Minimum value of stability correction
15  Real def_lambda_low = Real(0.001);
16  // Maximum value of stability correction
17  Real def_lambda_high = Real(0.04);
18  // Slope of change from lambda_low to lambda_high
19  Real def_lambda_slope = Real(2.65);
20  // stability threshold for which to apply more stability correction
21  Real def_lambda_thresh = Real(0.02);
22  // Temperature variance tuning factor
23  Real def_thl2tune = one;
24  // Moisture variance tuning factor
25  Real def_qw2tune = one;
26  // Temperature moisture covariance
27  Real def_qwthl2tune = one;
28  // Vertical velocity variance
29  Real def_w2tune = one;
30  // Length scale factor
31  Real def_length_fac = Real(0.5);
32  // Third moment vertical velocity damping factor
33  Real def_c_diag_3rd_mom = Real(7.0);
34  // Eddy diffusivity coefficient for heat
35  Real def_coeff_kh = Real(0.1);
36  // Eddy diffusivity coefficient for momentum
37  Real def_coeff_km = Real(0.1);
38 
39  runtime_options.lambda_low = def_lambda_low;
40  runtime_options.lambda_high = def_lambda_high;
41  runtime_options.lambda_slope = def_lambda_slope;
42  runtime_options.lambda_thresh = def_lambda_thresh;
43 
44  runtime_options.thl2tune = def_thl2tune;
45  runtime_options.qwthl2tune = def_qwthl2tune;
46  runtime_options.qw2tune = def_qw2tune;
47  runtime_options.w2tune = def_w2tune;
48 
49  runtime_options.length_fac = def_length_fac;
50  runtime_options.c_diag_3rd_mom = def_c_diag_3rd_mom;
51  runtime_options.Ckh = def_coeff_kh;
52  runtime_options.Ckm = def_coeff_km;
53  runtime_options.shoc_1p5tke = def_shoc_1p5tke;
54  runtime_options.extra_diags = extra_shoc_diags;
55 
56  // Construct parser object for following reads
57  ParmParse pp("erf.shoc");
58 
59  // Parse runtime inputs at start up
60  pp.query("lambda_low" , runtime_options.lambda_low );
61  pp.query("lambda_high" , runtime_options.lambda_high );
62  pp.query("lambda_slope" , runtime_options.lambda_slope );
63  pp.query("lambda_thresh" , runtime_options.lambda_thresh );
64  pp.query("thl2tune" , runtime_options.thl2tune );
65  pp.query("qw2tune" , runtime_options.qw2tune );
66  pp.query("qwthl2tune" , runtime_options.qwthl2tune );
67  pp.query("w2tune" , runtime_options.w2tune );
68  pp.query("length_fac" , runtime_options.length_fac );
69  pp.query("c_diag_3rd_mom" , runtime_options.c_diag_3rd_mom);
70  pp.query("coeff_kh" , runtime_options.Ckh );
71  pp.query("coeff_km" , runtime_options.Ckm );
72  pp.query("shoc_1p5tke" , runtime_options.shoc_1p5tke );
73  pp.query("extra_shoc_diags", runtime_options.extra_diags );
74 
75  // Set to default but allow us to change it through the inputs file
76  pp.query("apply_tms", apply_tms);
77  pp.query("check_flux_state", check_flux_state);
78  pp.query("extra_shoc_diags", extra_shoc_diags);
79  pp.query("column_conservation_check", column_conservation_check);
80 }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
ParmParse pp("prob")
amrex::Real Real
Definition: ERF_ShocInterface.H:19
bool check_flux_state
Definition: ERF_ShocInterface.H:648
bool extra_shoc_diags
Definition: ERF_ShocInterface.H:649
bool apply_tms
Definition: ERF_ShocInterface.H:647
SHF::SHOCRuntime runtime_options
Definition: ERF_ShocInterface.H:629
bool column_conservation_check
Definition: ERF_ShocInterface.H:650
Here is the call graph for this function:

Member Function Documentation

◆ add_fast_tend()

void SHOCInterface::add_fast_tend ( amrex::Vector< amrex::MultiFab > &  S_rhs)
704 {
705  for (MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
706  const auto& vbx_cc = mfi.validbox();
707  const auto& vbx_x = convert(vbx_cc,IntVect(1,0,0));
708  const auto& vbx_y = convert(vbx_cc,IntVect(0,1,0));
709 
710  const Array4<const Real>& c_arr = m_cons->const_array(mfi);
711 
712  const Array4<Real>& cc_rhs_arr = S_rhs[IntVars::cons].array(mfi);
713  const Array4<Real>& ru_rhs_arr = S_rhs[IntVars::xmom].array(mfi);
714  const Array4<Real>& rv_rhs_arr = S_rhs[IntVars::ymom].array(mfi);
715 
716  const Array4<const Real>& c_tend_arr = c_tend.const_array(mfi);
717  const Array4<const Real>& u_tend_arr = u_tend.const_array(mfi);
718  const Array4<const Real>& v_tend_arr = v_tend.const_array(mfi);
719 
720  ParallelFor(vbx_cc, vbx_x, vbx_y,
721  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
722  {
723  cc_rhs_arr(i,j,k,RhoTheta_comp) += c_arr(i,j,k,Rho_comp) * c_tend_arr(i,j,k,RhoTheta_comp);
724  },
725  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
726  {
727  ru_rhs_arr(i,j,k) += c_arr(i,j,k,Rho_comp) * u_tend_arr(i,j,k);
728  },
729  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
730  {
731  rv_rhs_arr(i,j,k) += c_arr(i,j,k,Rho_comp) * v_tend_arr(i,j,k);
732  });
733  }
734 }
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
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+myhalf) *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<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;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::MultiFab c_tend
Definition: ERF_ShocInterface.H:614
amrex::MultiFab * m_cons
Definition: ERF_ShocInterface.H:591
amrex::MultiFab u_tend
Definition: ERF_ShocInterface.H:615
amrex::MultiFab v_tend
Definition: ERF_ShocInterface.H:616
@ ymom
Definition: ERF_IndexDefines.H:178
@ cons
Definition: ERF_IndexDefines.H:176
@ xmom
Definition: ERF_IndexDefines.H:177
Here is the call graph for this function:

◆ add_slow_tend()

void SHOCInterface::add_slow_tend ( const amrex::MFIter &  mfi,
const amrex::Box &  tbx,
const amrex::Array4< amrex::Real > &  cc_rhs_arr 
)
741 {
742  bool moist = (m_cons->nComp() > RhoQ1_comp);
743 
744  const Array4<const Real>& c_arr = m_cons->const_array(mfi);
745 
746  const Array4<const Real>& c_tend_arr = c_tend.const_array(mfi);
747 
748  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
749  {
750  cc_rhs_arr(i,j,k,RhoKE_comp) += c_arr(i,j,k,Rho_comp) * c_tend_arr(i,j,k,RhoKE_comp);
751  if (moist) {
752  cc_rhs_arr(i,j,k,RhoQ1_comp) += c_arr(i,j,k,Rho_comp) * c_tend_arr(i,j,k,RhoQ1_comp);
753  cc_rhs_arr(i,j,k,RhoQ2_comp) += c_arr(i,j,k,Rho_comp) * c_tend_arr(i,j,k,RhoQ2_comp);
754  }
755  });
756 }
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
Here is the call graph for this function:

◆ alloc_buffers()

void SHOCInterface::alloc_buffers ( )
200 {
201  // Interface data structures
202  //=======================================================
203  omega = view_2d("Omega" , m_num_cols, m_num_layers );
204  surf_sens_flux = view_1d("Sfc sens flux" , m_num_cols);
205  surf_mom_flux = sview_2d("Sfc mom flux" , m_num_cols, m_num_vel_comp);
206  surf_evap = view_1d("Sfc evap" , m_num_cols);
207  T_mid = view_2d("T_mid" , m_num_cols, m_num_layers );
208  qv = view_2d("Qv" , m_num_cols, m_num_layers );
209  surf_drag_coeff_tms = view_1d("surf_drag_coeff", m_num_cols);
210 
211  // Input data structures
212  //=======================================================
213  p_mid = view_2d("P_mid" , m_num_cols, m_num_layers );
214  p_int = view_2d("P_int" , m_num_cols, m_num_layers+1);
215  pseudo_dens = view_2d("Pseudo density" , m_num_cols, m_num_layers );
216  phis = view_1d("Phis" , m_num_cols);
217 
218  // Input/Output data structures
219  //=======================================================
221  cldfrac_liq = view_2d("Cld_frac_liq" , m_num_cols, m_num_layers );
222  tke = view_2d("Tke" , m_num_cols, m_num_layers );
223  qc = view_2d("Qc" , m_num_cols, m_num_layers );
224 
225  // Output data structures
226  //=======================================================
227  pblh = view_1d("pbl_height" , m_num_cols);
228  inv_qc_relvar = view_2d("inv_qc_relvar" , m_num_cols, m_num_layers);
229  tkh = view_2d("eddy_diff_heat" , m_num_cols, m_num_layers);
230  w_sec = view_2d("w_sec" , m_num_cols, m_num_layers);
231  cldfrac_liq_prev = view_2d("cld_frac_liq_prev" , m_num_cols, m_num_layers );
232  ustar = view_1d("ustar" , m_num_cols);
233  obklen = view_1d("obklen" , m_num_cols);
234 
235  // Extra diagnostic data structures
236  //=======================================================
237  if (extra_shoc_diags) {
238  brunt = view_2d("brunt" , m_num_cols, m_num_layers);
239  shoc_mix = view_2d("shoc_mix" , m_num_cols, m_num_layers);
240  isotropy = view_2d("isotropy" , m_num_cols, m_num_layers);
241  shoc_cond = view_2d("shoc_cond", m_num_cols, m_num_layers);
242  shoc_evap = view_2d("shoc_evap", m_num_cols, m_num_layers);
243 
244  wthl_sec = view_2d("wthl_sec" , m_num_cols, m_num_layers+1);
245  thl_sec = view_2d("thl_sec" , m_num_cols, m_num_layers+1);
246  wqw_sec = view_2d("wqw_sec" , m_num_cols, m_num_layers+1);
247  qw_sec = view_2d("qw_sec" , m_num_cols, m_num_layers+1);
248  uw_sec = view_2d("uw_sec" , m_num_cols, m_num_layers+1);
249  vw_sec = view_2d("vw_sec" , m_num_cols, m_num_layers+1);
250  w3 = view_2d("w3" , m_num_cols, m_num_layers+1);
251  }
252 
253  // Tracer data structures
254  //=======================================================
255  // NOTE: Use layoutright format
256  Kokkos::LayoutStride layout(m_num_cols , m_num_cols*m_num_layers, // stride for dim0
257  m_num_layers , m_num_layers, // stride for dim1
258  m_num_tracers, 1 // stride for dim2
259  );
260  qtracers = view_3d_strided("Qtracers" , layout);
261 
262  // Boundary flux data structures
263  //=======================================================
265  vapor_flux = view_1d("vapor_flux", m_num_cols);
266  water_flux = view_1d("water_flux", m_num_cols);
267  ice_flux = view_1d("ice_flux" , m_num_cols);
268  heat_flux = view_1d("heat_flux" , m_num_cols);
269  }
270 }
view_2d p_int
Definition: ERF_ShocInterface.H:665
view_2d wqw_sec
Definition: ERF_ShocInterface.H:698
view_1d surf_drag_coeff_tms
Definition: ERF_ShocInterface.H:660
view_2d thl_sec
Definition: ERF_ShocInterface.H:697
view_2d tkh
Definition: ERF_ShocInterface.H:682
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:674
typename KT::template view_2d< Real > sview_2d
Definition: ERF_ShocInterface.H:51
Int m_num_layers
Definition: ERF_ShocInterface.H:571
view_2d vw_sec
Definition: ERF_ShocInterface.H:701
view_1d ice_flux
Definition: ERF_ShocInterface.H:712
view_2d shoc_cond
Definition: ERF_ShocInterface.H:693
view_2d brunt
Definition: ERF_ShocInterface.H:690
view_1d heat_flux
Definition: ERF_ShocInterface.H:713
Int m_num_vel_comp
Definition: ERF_ShocInterface.H:575
view_1d surf_sens_flux
Definition: ERF_ShocInterface.H:655
view_1d pblh
Definition: ERF_ShocInterface.H:680
Int m_num_cols
Definition: ERF_ShocInterface.H:570
view_2d p_mid
Definition: ERF_ShocInterface.H:664
view_2d uw_sec
Definition: ERF_ShocInterface.H:700
view_2d qc
Definition: ERF_ShocInterface.H:676
view_2d shoc_mix
Definition: ERF_ShocInterface.H:691
view_2d qw_sec
Definition: ERF_ShocInterface.H:699
view_2d cldfrac_liq_prev
Definition: ERF_ShocInterface.H:684
view_2d omega
Definition: ERF_ShocInterface.H:654
typename SHF::view_3d_strided< Spack > view_3d_strided
Definition: ERF_ShocInterface.H:55
view_1d vapor_flux
Definition: ERF_ShocInterface.H:710
view_1d obklen
Definition: ERF_ShocInterface.H:686
view_2d pseudo_dens
Definition: ERF_ShocInterface.H:666
view_2d inv_qc_relvar
Definition: ERF_ShocInterface.H:681
view_1d surf_evap
Definition: ERF_ShocInterface.H:657
view_2d w_sec
Definition: ERF_ShocInterface.H:683
view_1d phis
Definition: ERF_ShocInterface.H:667
view_2d shoc_evap
Definition: ERF_ShocInterface.H:694
view_3d horiz_wind
Definition: ERF_ShocInterface.H:671
typename SHF::view_1d< Real > view_1d
Definition: ERF_ShocInterface.H:47
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:706
view_1d ustar
Definition: ERF_ShocInterface.H:685
view_2d w3
Definition: ERF_ShocInterface.H:702
Int m_num_tracers
Definition: ERF_ShocInterface.H:574
sview_2d surf_mom_flux
Definition: ERF_ShocInterface.H:656
typename SHF::view_2d< SHF::Spack > view_2d
Definition: ERF_ShocInterface.H:49
view_2d qv
Definition: ERF_ShocInterface.H:659
view_2d wthl_sec
Definition: ERF_ShocInterface.H:696
view_2d tke
Definition: ERF_ShocInterface.H:675
view_1d water_flux
Definition: ERF_ShocInterface.H:711
view_2d T_mid
Definition: ERF_ShocInterface.H:658
typename SHF::view_3d< Spack > view_3d
Definition: ERF_ShocInterface.H:53
view_2d isotropy
Definition: ERF_ShocInterface.H:692

◆ apply_turbulent_mountain_stress()

void SHOCInterface::apply_turbulent_mountain_stress ( )
protected
1105 {
1106  auto rrho_i = m_buffer.rrho_i;
1107  auto upwp_sfc = m_buffer.upwp_sfc;
1108  auto vpwp_sfc = m_buffer.vpwp_sfc;
1109  auto surf_drag_coeff_tms_d = surf_drag_coeff_tms;
1110  auto horiz_wind_d = horiz_wind;
1111 
1112  const int nlev_v = (m_num_layers-1)/Spack::n;
1113  const int nlev_p = (m_num_layers-1)%Spack::n;
1114  const int nlevi_v = m_num_layers/Spack::n;
1115  const int nlevi_p = m_num_layers%Spack::n;
1116 
1117  Kokkos::parallel_for("apply_tms", KT::RangePolicy(0, m_num_cols), KOKKOS_LAMBDA (const int i)
1118  {
1119  upwp_sfc(i) -= surf_drag_coeff_tms_d(i)*horiz_wind_d(i,0,nlev_v)[nlev_p]/rrho_i(i,nlevi_v)[nlevi_p];
1120  vpwp_sfc(i) -= surf_drag_coeff_tms_d(i)*horiz_wind_d(i,1,nlev_v)[nlev_p]/rrho_i(i,nlevi_v)[nlevi_p];
1121  });
1122 }
Buffer m_buffer
Definition: ERF_ShocInterface.H:622
uview_1d< Real > upwp_sfc
Definition: ERF_ShocInterface.H:476
uview_1d< Real > vpwp_sfc
Definition: ERF_ShocInterface.H:477
uview_2d< Spack > rrho_i
Definition: ERF_ShocInterface.H:497

◆ check_flux_state_consistency()

void SHOCInterface::check_flux_state_consistency ( const double  dt)
protected
1127 {
1128  using PC = scream::physics::Constants<Real>;
1129  using RU = ekat::ReductionUtils<KT::ExeSpace>;
1130  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
1131 
1132  const Real gravit = PC::gravit;
1133  const Real qmin = 1e-12; // minimum permitted constituent concentration (kg/kg)
1134 
1135  const auto& pseudo_density = pseudo_dens;
1136  auto qv_d = qv;
1137  auto surf_evap_d = surf_evap;
1138 
1139  const auto nlevs = m_num_layers;
1140  const auto nlev_packs = ekat::npack<Spack>(nlevs);
1141  const auto last_pack_idx = (nlevs-1)/Spack::n;
1142  const auto last_pack_entry = (nlevs-1)%Spack::n;
1143  const auto policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
1144  Kokkos::parallel_for("check_flux_state_consistency",
1145  policy,
1146  KOKKOS_LAMBDA (const KT::MemberType& team)
1147  {
1148  const auto i = team.league_rank();
1149 
1150  const auto& pseudo_density_i = ekat::subview(pseudo_density, i);
1151  const auto& qv_i = ekat::subview(qv_d, i);
1152 
1153  // reciprocal of pseudo_density at the bottom layer
1154  const auto rpdel = one/pseudo_density_i(last_pack_idx)[last_pack_entry];
1155 
1156  // Check if the negative surface latent heat flux can exhaust
1157  // the moisture in the lowest model level. If so, apply fixer.
1158  const auto condition = surf_evap_d(i) - (qmin - qv_i(last_pack_idx)[last_pack_entry])/(dt*gravit*rpdel);
1159  if (condition < 0) {
1160  const auto cc = std::fabs(surf_evap_d(i)*dt*gravit);
1161 
1162  auto tracer_mass = [&](const int k)
1163  {
1164  return qv_i(k)*pseudo_density_i(k);
1165  };
1166  Real mm = RU::view_reduction(team, 0, nlevs, tracer_mass);
1167 
1168  EKAT_KERNEL_ASSERT_MSG(mm >= cc, "Error! Total mass of column vapor should be greater than mass of surf_evap.\n");
1169 
1170  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&](const int& k)
1171  {
1172  const auto adjust = cc*qv_i(k)*pseudo_density_i(k)/mm;
1173  qv_i(k) = (qv_i(k)*pseudo_density_i(k) - adjust)/pseudo_density_i(k);
1174  });
1175 
1176  surf_evap_d(i) = 0;
1177  }
1178  });
1179 }
@ PC
Definition: ERF_IndexDefines.H:123

◆ dealloc_buffers()

void SHOCInterface::dealloc_buffers ( )
275 {
276  // Contiguous memory buffer view
277  //=======================================================
279 
280  // Interface data structures
281  //=======================================================
282  omega = view_2d();
285  surf_evap = view_1d();
286  T_mid = view_2d();
287  qv = view_2d();
289 
290  // Input data structures
291  //=======================================================
292  p_mid = view_2d();
293  p_int = view_2d();
294  pseudo_dens = view_2d();
295  phis = view_1d();
296 
297  // Input/Output data structures
298  //=======================================================
299  horiz_wind = view_3d();
300  cldfrac_liq = view_2d();
301  tke = view_2d();
302  qc = view_2d();
303 
304  // Output data structures
305  //=======================================================
306  pblh = view_1d();
308  tkh = view_2d();
309  w_sec = view_2d();
311  ustar = view_1d();
312  obklen = view_1d();
313 
314  // Extra diagnostic data structures
315  //=======================================================
316  if (extra_shoc_diags) {
317  brunt = view_2d();
318  shoc_mix = view_2d();
319  isotropy = view_2d();
320  shoc_cond = view_2d();
321  shoc_evap = view_2d();
322 
323  wthl_sec = view_2d();
324  thl_sec = view_2d();
325  wqw_sec = view_2d();
326  qw_sec = view_2d();
327  uw_sec = view_2d();
328  vw_sec = view_2d();
329  w3 = view_2d();
330  }
331 
332  // Tracer data structures
333  //=======================================================
335 
336  // Boundary flux data structures
337  //=======================================================
339  vapor_flux = view_1d();
340  water_flux = view_1d();
341  ice_flux = view_1d();
342  heat_flux = view_1d();
343  }
344 }
view_1d tot_buff_view
Definition: ERF_ShocInterface.H:643

◆ finalize_impl()

void SHOCInterface::finalize_impl ( const Real  dt)
1092 {
1093  // Do nothing (per SHOCMacrophysics::finalize_impl())
1094 
1095  // Fill the AMReX MFs from Kokkos Views
1097 
1098  // Deallocate the buffer arrays
1099  dealloc_buffers();
1100 }
void kokkos_buffers_to_mf(const Real dt)
Definition: ERF_ShocInterface.cpp:517
void dealloc_buffers()
Definition: ERF_ShocInterface.cpp:274

◆ init_buffers()

void SHOCInterface::init_buffers ( )
protected
787 {
788 
789  // Buffer of contiguous memory
790  auto buffer_size = requested_buffer_size_in_bytes();
791  tot_buff_view = view_1d("contiguous shoc_buffer",buffer_size);
792  Real* mem = reinterpret_cast<Real*>(tot_buff_view.data());
793 
794  // 1d scalar views
795  using scalar_view_t = decltype(m_buffer.wpthlp_sfc);
796  scalar_view_t* _1d_scalar_view_ptrs[Buffer::num_1d_scalar_ncol] =
798 #ifdef SCREAM_SHOC_SMALL_KERNELS
799  , &m_buffer.se_b, &m_buffer.ke_b, &m_buffer.wv_b, &m_buffer.wl_b
800  , &m_buffer.se_a, &m_buffer.ke_a, &m_buffer.wv_a, &m_buffer.wl_a
801  , &m_buffer.kbfs, &m_buffer.ustar2, &m_buffer.wstar
802 #endif
803  };
804  for (int i = 0; i < Buffer::num_1d_scalar_ncol; ++i) {
805  *_1d_scalar_view_ptrs[i] = scalar_view_t(mem, m_num_cols);
806  mem += _1d_scalar_view_ptrs[i]->size();
807  }
808 
809  Spack* s_mem = reinterpret_cast<Spack*>(mem);
810 
811  // 2d packed views
812  const int nlev_packs = ekat::npack<Spack>(m_num_layers);
813  const int nlevi_packs = ekat::npack<Spack>(m_num_layers+1);
814  const int num_tracer_packs = ekat::npack<Spack>(m_num_tracers);
815 
816  m_buffer.pref_mid = decltype(m_buffer.pref_mid)(s_mem, nlev_packs);
817  s_mem += m_buffer.pref_mid.size();
818 
819  using spack_2d_view_t = decltype(m_buffer.z_mid);
820  spack_2d_view_t* _2d_spack_mid_view_ptrs[Buffer::num_2d_vector_mid] =
824 #ifdef SCREAM_SHOC_SMALL_KERNELS
825  , &m_buffer.rho_zt, &m_buffer.shoc_qv, &m_buffer.tabs, &m_buffer.dz_zt
826 #endif
827  };
828 
829  spack_2d_view_t* _2d_spack_int_view_ptrs[Buffer::num_2d_vector_int] =
833 #ifdef SCREAM_SHOC_SMALL_KERNELS
834  , &m_buffer.dz_zi
835 #endif
836  };
837 
838  for (int i = 0; i < Buffer::num_2d_vector_mid; ++i) {
839  *_2d_spack_mid_view_ptrs[i] = spack_2d_view_t(s_mem, m_num_cols, nlev_packs);
840  s_mem += _2d_spack_mid_view_ptrs[i]->size();
841  }
842 
843  for (int i = 0; i < Buffer::num_2d_vector_int; ++i) {
844  *_2d_spack_int_view_ptrs[i] = spack_2d_view_t(s_mem, m_num_cols, nlevi_packs);
845  s_mem += _2d_spack_int_view_ptrs[i]->size();
846  }
847  m_buffer.wtracer_sfc = decltype(m_buffer.wtracer_sfc)(s_mem, m_num_cols, num_tracer_packs);
848  s_mem += m_buffer.wtracer_sfc.size();
849 
850  // WSM data
851  m_buffer.wsm_data = s_mem;
852 }
size_t requested_buffer_size_in_bytes() const
Definition: ERF_ShocInterface.cpp:760
typename SHF::Spack Spack
Definition: ERF_ShocInterface.H:43
uview_2d< Spack > w3
Definition: ERF_ShocInterface.H:522
uview_2d< Spack > qc_copy
Definition: ERF_ShocInterface.H:509
uview_2d< Spack > dse
Definition: ERF_ShocInterface.H:507
uview_2d< Spack > vw_sec
Definition: ERF_ShocInterface.H:521
uview_1d< Real > wpthlp_sfc
Definition: ERF_ShocInterface.H:474
uview_2d< Spack > thv
Definition: ERF_ShocInterface.H:498
uview_2d< Spack > wqls_sec
Definition: ERF_ShocInterface.H:523
uview_2d< Spack > wthl_sec
Definition: ERF_ShocInterface.H:517
uview_2d< Spack > thlm
Definition: ERF_ShocInterface.H:505
uview_2d< Spack > uw_sec
Definition: ERF_ShocInterface.H:520
uview_2d< Spack > z_mid
Definition: ERF_ShocInterface.H:494
uview_2d< Spack > rrho
Definition: ERF_ShocInterface.H:496
uview_2d< Spack > wtke_sec
Definition: ERF_ShocInterface.H:519
uview_2d< Spack > unused
Definition: ERF_ShocInterface.H:493
Spack * wsm_data
Definition: ERF_ShocInterface.H:534
static constexpr int num_2d_vector_int
Definition: ERF_ShocInterface.H:467
uview_2d< Spack > wtracer_sfc
Definition: ERF_ShocInterface.H:502
uview_2d< Spack > inv_exner
Definition: ERF_ShocInterface.H:504
uview_2d< Spack > qwthl_sec
Definition: ERF_ShocInterface.H:516
uview_2d< Spack > qw
Definition: ERF_ShocInterface.H:506
uview_2d< Spack > brunt
Definition: ERF_ShocInterface.H:524
uview_2d< Spack > shoc_mix
Definition: ERF_ShocInterface.H:511
uview_2d< Spack > tke_copy
Definition: ERF_ShocInterface.H:508
static constexpr int num_1d_scalar_ncol
Definition: ERF_ShocInterface.H:460
uview_2d< Spack > isotropy
Definition: ERF_ShocInterface.H:512
uview_2d< Spack > zt_grid
Definition: ERF_ShocInterface.H:500
uview_2d< Spack > w_sec
Definition: ERF_ShocInterface.H:513
uview_2d< Spack > thl_sec
Definition: ERF_ShocInterface.H:514
uview_2d< Spack > wm_zt
Definition: ERF_ShocInterface.H:503
uview_2d< Spack > qw_sec
Definition: ERF_ShocInterface.H:515
uview_2d< Spack > wqw_sec
Definition: ERF_ShocInterface.H:518
uview_2d< Spack > z_int
Definition: ERF_ShocInterface.H:495
uview_2d< Spack > shoc_ql2
Definition: ERF_ShocInterface.H:510
uview_1d< Real > wprtp_sfc
Definition: ERF_ShocInterface.H:475
uview_2d< Spack > zi_grid
Definition: ERF_ShocInterface.H:501
static constexpr int num_2d_vector_mid
Definition: ERF_ShocInterface.H:466
uview_1d< Spack > pref_mid
Definition: ERF_ShocInterface.H:491
uview_2d< Spack > dz
Definition: ERF_ShocInterface.H:499

◆ initialize_impl()

void SHOCInterface::initialize_impl ( )
857 {
858  // Alias local variables from temporary buffer
859  auto z_mid = m_buffer.z_mid;
860  auto z_int = m_buffer.z_int;
861  auto wpthlp_sfc = m_buffer.wpthlp_sfc;
862  auto wprtp_sfc = m_buffer.wprtp_sfc;
863  auto upwp_sfc = m_buffer.upwp_sfc;
864  auto vpwp_sfc = m_buffer.vpwp_sfc;
865  auto rrho = m_buffer.rrho;
866  auto rrho_i = m_buffer.rrho_i;
867  auto thv = m_buffer.thv;
868  auto dz = m_buffer.dz;
869  auto zt_grid = m_buffer.zt_grid;
870  auto zi_grid = m_buffer.zi_grid;
871  auto wtracer_sfc = m_buffer.wtracer_sfc;
872  auto wm_zt = m_buffer.wm_zt;
873  auto inv_exner = m_buffer.inv_exner;
874  auto thlm = m_buffer.thlm;
875  auto qw = m_buffer.qw;
876  auto dse = m_buffer.dse;
877  auto tke_copy = m_buffer.tke_copy;
878  auto qc_copy = m_buffer.qc_copy;
879  auto shoc_ql2 = m_buffer.shoc_ql2;
880 
881  // For now, set z_int(i,nlevs) = z_surf = 0
882  const Real z_surf = Real(0.);
883 
884  // Set preprocess variables
887  surf_mom_flux, qtracers, qv, qc, qc_copy, tke, tke_copy, z_mid, z_int,
888  dse, rrho, rrho_i, thv, dz, zt_grid, zi_grid, wpthlp_sfc, wprtp_sfc, upwp_sfc,
889  vpwp_sfc, wtracer_sfc, wm_zt, inv_exner, thlm, qw, cldfrac_liq, cldfrac_liq_prev);
890 
891  // Input Variables:
892  input.zt_grid = shoc_preprocess.zt_grid;
893  input.zi_grid = shoc_preprocess.zi_grid;
894  input.pres = p_mid;
895  input.presi = p_int;
896  input.pdel = pseudo_dens;
897  input.thv = shoc_preprocess.thv;
898  input.w_field = shoc_preprocess.wm_zt;
899  input.wthl_sfc = shoc_preprocess.wpthlp_sfc;
900  input.wqw_sfc = shoc_preprocess.wprtp_sfc;
901  input.uw_sfc = shoc_preprocess.upwp_sfc;
902  input.vw_sfc = shoc_preprocess.vpwp_sfc;
903  input.wtracer_sfc = shoc_preprocess.wtracer_sfc;
904  input.inv_exner = shoc_preprocess.inv_exner;
905  input.phis = phis;
906 
907  // Input/Output Variables
912  input_output.horiz_wind = horiz_wind;
913  input_output.wthv_sec = sgs_buoy_flux;
915  input_output.tk = tk;
916  input_output.shoc_cldfrac = cldfrac_liq;
917  input_output.shoc_ql = qc_copy;
918 
919  // Output Variables
920  output.pblh = pblh;
921  output.shoc_ql2 = shoc_ql2;
922  output.tkh = tkh;
923  output.ustar = ustar;
924  output.obklen = obklen;
925 
926  // Output (diagnostic)
927  history_output.shoc_mix = m_buffer.shoc_mix;
928  history_output.isotropy = m_buffer.isotropy;
929  if (extra_shoc_diags) {
930  history_output.shoc_cond = shoc_cond;
931  history_output.shoc_evap = shoc_evap;
932  } else {
933  history_output.shoc_cond = m_buffer.unused;
934  history_output.shoc_evap = m_buffer.unused;
935  }
936  history_output.w_sec = w_sec;
937  history_output.thl_sec = m_buffer.thl_sec;
938  history_output.qw_sec = m_buffer.qw_sec;
939  history_output.qwthl_sec = m_buffer.qwthl_sec;
940  history_output.wthl_sec = m_buffer.wthl_sec;
941  history_output.wqw_sec = m_buffer.wqw_sec;
942  history_output.wtke_sec = m_buffer.wtke_sec;
943  history_output.uw_sec = m_buffer.uw_sec;
944  history_output.vw_sec = m_buffer.vw_sec;
946  history_output.wqls_sec = m_buffer.wqls_sec;
947  history_output.brunt = m_buffer.brunt;
948 
949 #ifdef SCREAM_SHOC_SMALL_KERNELS
950  temporaries.se_b = m_buffer.se_b;
951  temporaries.ke_b = m_buffer.ke_b;
952  temporaries.wv_b = m_buffer.wv_b;
953  temporaries.wl_b = m_buffer.wl_b;
954  temporaries.se_a = m_buffer.se_a;
955  temporaries.ke_a = m_buffer.ke_a;
956  temporaries.wv_a = m_buffer.wv_a;
957  temporaries.wl_a = m_buffer.wl_a;
958  temporaries.kbfs = m_buffer.kbfs;
959  temporaries.ustar2 = m_buffer.ustar2;
960  temporaries.wstar = m_buffer.wstar;
961 
962  temporaries.rho_zt = m_buffer.rho_zt;
963  temporaries.shoc_qv = m_buffer.shoc_qv;
964  temporaries.tabs = m_buffer.tabs;
965  temporaries.dz_zt = m_buffer.dz_zt;
966  temporaries.dz_zi = m_buffer.dz_zi;
967 #endif
968 
969  // Set postprocess variables
971  rrho, qv, qw, qc, qc_copy, tke,
972  tke_copy, qtracers, shoc_ql2,
974  T_mid, dse, z_mid, phis);
975 
976  // Setup WSM for internal local variables
977  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
978  const auto nlev_packs = ekat::npack<Spack>(m_num_layers);
979  const auto nlevi_packs = ekat::npack<Spack>(m_num_layers+1);
980  const int n_wind_slots = ekat::npack<Spack>(m_num_vel_comp)*Spack::n;
981  const int n_trac_slots = ekat::npack<Spack>(m_num_tracers)*Spack::n;
982  const auto default_policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
983  workspace_mgr.setup(m_buffer.wsm_data, nlevi_packs, 14+(n_wind_slots+n_trac_slots), default_policy);
984 
985  // NOTE: Vertical indices were permuted, so top and bottom are correct
986  // Maximum number of levels in pbl from surface
987  const int ntop_shoc = 0;
988  const int nbot_shoc = m_num_layers;
989  auto p_mid_d = p_mid;
990  view_1d pref_mid("pref_mid", m_num_layers);
991  Spack* s_mem = reinterpret_cast<Spack*>(pref_mid.data());
992  SHF::view_1d<Spack> pref_mid_um(s_mem, m_num_layers);
993  const auto policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
994  Kokkos::parallel_for("pref_mid",
995  policy,
996  KOKKOS_LAMBDA (const KT::MemberType& team)
997  {
998  const auto i = team.league_rank();
999  if (i==0) {
1000  const auto& pmid_i = ekat::subview(p_mid_d, i);
1001  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&](const int& k)
1002  {
1003  pref_mid_um(k) = pmid_i(k);
1004  });
1005  }
1006  });
1007  Kokkos::fence();
1008  m_npbl = SHF::shoc_init(nbot_shoc,ntop_shoc,pref_mid_um);
1009 
1010  // Cell length for input dx and dy
1011  view_1d cell_length_x("cell_length_x", m_num_cols);
1012  view_1d cell_length_y("cell_length_y", m_num_cols);
1013  Kokkos::deep_copy(cell_length_x, m_geom.CellSize(0));
1014  Kokkos::deep_copy(cell_length_y, m_geom.CellSize(1));
1015  input.dx = cell_length_x;
1016  input.dy = cell_length_y;
1017 }
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
SHOCPreprocess shoc_preprocess
Definition: ERF_ShocInterface.H:635
amrex::Geometry m_geom
Definition: ERF_ShocInterface.H:585
view_2d sgs_buoy_flux
Definition: ERF_ShocInterface.H:672
SHF::SHOCOutput output
Definition: ERF_ShocInterface.H:627
ekat::WorkspaceManager< Spack, KT::Device > workspace_mgr
Definition: ERF_ShocInterface.H:639
Int m_npbl
Definition: ERF_ShocInterface.H:572
SHF::SHOCInput input
Definition: ERF_ShocInterface.H:625
SHOCPostprocess shoc_postprocess
Definition: ERF_ShocInterface.H:636
SHF::SHOCInputOutput input_output
Definition: ERF_ShocInterface.H:626
view_2d tk
Definition: ERF_ShocInterface.H:673
SHF::SHOCHistoryOutput history_output
Definition: ERF_ShocInterface.H:628
void set_variables(const int ncol_, const int nlev_, const view_2d_const &rrho_, const view_2d &qv_, const view_2d_const &qw_, const view_2d &qc_, const view_2d_const &qc_copy_, const view_2d &tke_, const view_2d_const &tke_copy_, const view_3d_strided &qtracers_, const view_2d_const &qc2_, const view_2d &cldfrac_liq_, const view_2d &inv_qc_relvar_, const view_2d &T_mid_, const view_2d_const &dse_, const view_2d_const &z_mid_, const view_1d_const phis_)
Definition: ERF_ShocInterface.H:416
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:252
view_2d thlm
Definition: ERF_ShocInterface.H:274
view_2d zt_grid
Definition: ERF_ShocInterface.H:265
view_2d inv_exner
Definition: ERF_ShocInterface.H:273
view_2d wm_zt
Definition: ERF_ShocInterface.H:272
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:271
view_1d wprtp_sfc
Definition: ERF_ShocInterface.H:268
view_2d qw
Definition: ERF_ShocInterface.H:275
view_1d upwp_sfc
Definition: ERF_ShocInterface.H:269
view_1d wpthlp_sfc
Definition: ERF_ShocInterface.H:267
view_1d vpwp_sfc
Definition: ERF_ShocInterface.H:270
view_2d shoc_s
Definition: ERF_ShocInterface.H:258
view_2d tke_copy
Definition: ERF_ShocInterface.H:260
view_2d zi_grid
Definition: ERF_ShocInterface.H:266
view_2d thv
Definition: ERF_ShocInterface.H:263
void set_variables(const int ncol_, const int nlev_, const Real z_surf_, const view_2d_const &T_mid_, const view_2d_const &p_mid_, const view_2d_const &p_int_, const view_2d_const &pseudo_density_, const view_2d_const &omega_, const view_1d_const &phis_, const view_1d_const &surf_sens_flux_, const view_1d_const &surf_evap_, const sview_2d_const &surf_mom_flux_, const view_3d_strided &qtracers_, const view_2d &qv_, const view_2d_const &qc_, const view_2d &qc_copy_, const view_2d &tke_, const view_2d &tke_copy_, const view_2d &z_mid_, const view_2d &z_int_, const view_2d &dse_, const view_2d &rrho_, const view_2d &rrho_i_, const view_2d &thv_, const view_2d &dz_, const view_2d &zt_grid_, const view_2d &zi_grid_, const view_1d &wpthlp_sfc_, const view_1d &wprtp_sfc_, const view_1d &upwp_sfc_, const view_1d &vpwp_sfc_, const view_2d &wtracer_sfc_, const view_2d &wm_zt_, const view_2d &inv_exner_, const view_2d &thlm_, const view_2d &qw_, const view_2d &cldfrac_liq_, const view_2d &cldfrac_liq_prev_)
Definition: ERF_ShocInterface.H:281

◆ kokkos_buffers_to_mf()

void SHOCInterface::kokkos_buffers_to_mf ( const Real  dt)
518 {
519  //
520  // Expose for device capture
521  //
522 
523  // Buffer data structures
524  //=======================================================
525  auto thlm_d = m_buffer.thlm;
526 
527  // Interface data structures
528  //=======================================================
529  auto T_mid_d = T_mid;
530  auto qv_d = qv;
531 
532  // Input/Output data structures
533  //=======================================================
534  auto horiz_wind_d = horiz_wind;
535  auto tke_d = tke;
536  auto qc_d = qc;
537 
538  bool moist = (m_cons->nComp() > RhoQ1_comp);
539  for (MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
540  // NOTE: No ghost cells when going back to MFs
541  const auto& vbx_cc = mfi.validbox();
542  const auto& vbx_x = convert(vbx_cc,IntVect(1,0,0));
543  const auto& vbx_y = convert(vbx_cc,IntVect(0,1,0));
544 
545  // NOTE: Grown box only for mapping
546  const auto& gbx = mfi.tilebox(IntVect(0,0,0),IntVect(1,1,0));
547  const int nx = gbx.length(0);
548  const int imin = gbx.smallEnd(0);
549  const int jmin = gbx.smallEnd(1);
550  const int kmax = gbx.bigEnd(2);
551  const int offset = m_col_offsets[mfi.index()];
552 
553  const Array4<const Real>& cons_arr = m_cons->const_array(mfi);
554  const Array4<const Real>& u_arr = m_xvel->const_array(mfi);
555  const Array4<const Real>& v_arr = m_yvel->const_array(mfi);
556 
557  const Array4<Real>& c_tend_arr = c_tend.array(mfi);
558  const Array4<Real>& u_tend_arr = u_tend.array(mfi);
559  const Array4<Real>& v_tend_arr = v_tend.array(mfi);
560 
561  ParallelFor(vbx_cc, vbx_x, vbx_y,
562  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
563  {
564  // NOTE: k gets permuted with ilay
565  // map [i,j,k] 0-based to [icol, ilay] 0-based
566  const int icol = (j-jmin)*nx + (i-imin) + offset;
567  const int ilay = kmax - k;
568 
569  // Density at CC
570  Real r = cons_arr(i,j,k,Rho_comp);
571 
572  // Theta at CC (eamxx_common_physics_functions_impl.hpp L123)
573  Real Th = thlm_d(icol,ilay)[0] / ( one - (one / T_mid_d(icol,ilay)[0])
574  * (C::LatVap/C::Cpair) * qc_d(icol,ilay)[0] );
575 
576  // Populate the tendencies
577  c_tend_arr(i,j,k,RhoTheta_comp) = ( Th - cons_arr(i,j,k,RhoTheta_comp)/r ) / dt;
578  c_tend_arr(i,j,k,RhoKE_comp) = ( tke_d(icol,ilay)[0] - cons_arr(i,j,k,RhoKE_comp )/r ) / dt;
579  if (moist) {
580  c_tend_arr(i,j,k,RhoQ1_comp) = ( qv_d(icol,ilay)[0] - cons_arr(i,j,k,RhoQ1_comp)/r ) / dt;
581  c_tend_arr(i,j,k,RhoQ2_comp) = ( qc_d(icol,ilay)[0] - cons_arr(i,j,k,RhoQ2_comp)/r ) / dt;
582  }
583  },
584  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
585  {
586  // NOTE: k gets permuted with ilay
587  // map [i,j,k] 0-based to [icol, ilay] 0-based
588  const int icol = (j-jmin)*nx + (i-imin) + offset;
589  const int ilay = kmax - k;
590 
591  int icolim = (j-jmin)*nx + (i-1-imin) + offset;
592  Real uvel = Real(0.5) * (horiz_wind_d(icol,0,ilay)[0] + horiz_wind_d(icolim,0,ilay)[0]);
593  u_tend_arr(i,j,k) = ( uvel - u_arr(i,j,k) ) / dt;
594  },
595  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
596  {
597  // NOTE: k gets permuted with ilay
598  // map [i,j,k] 0-based to [icol, ilay] 0-based
599  const int icol = (j-jmin)*nx + (i-imin) + offset;
600  const int ilay = kmax - k;
601 
602  int icoljm = (j-1-jmin)*nx + (i-imin) + offset;
603  Real vvel = Real(0.5) * (horiz_wind_d(icol,1,ilay)[0] + horiz_wind_d(icoljm,1,ilay)[0]);
604  v_tend_arr(i,j,k) = ( vvel - v_arr(i,j,k) ) / dt;
605  });
606  }
607 }
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::MultiFab * m_yvel
Definition: ERF_ShocInterface.H:595
amrex::MultiFab * m_xvel
Definition: ERF_ShocInterface.H:594
amrex::Vector< int > m_col_offsets
Definition: ERF_ShocInterface.H:566
Here is the call graph for this function:

◆ mf_to_kokkos_buffers()

void SHOCInterface::mf_to_kokkos_buffers ( )
349 {
350  // FillBoundary for internal ghost cells for u/v averaging
351  m_tau13->FillBoundary(m_geom.periodicity());
352  m_tau23->FillBoundary(m_geom.periodicity());
353  m_hfx3->FillBoundary(m_geom.periodicity());
354  m_qfx3->FillBoundary(m_geom.periodicity());
355 
356  //
357  // Expose for device capture
358  //
359 
360  // Interface data structures
361  //=======================================================
362  auto omega_d = omega;
363  auto surf_sens_flux_d = surf_sens_flux;
364  auto surf_mom_flux_d = surf_mom_flux;
365  auto surf_evap_d = surf_evap;
366  auto T_mid_d = T_mid;
367  auto qv_d = qv;
368  auto surf_drag_coeff_tms_d = surf_drag_coeff_tms;
369 
370  // Input data structures
371  //=======================================================
372  auto p_mid_d = p_mid;
373  auto p_int_d = p_int;
374  auto pseudo_dens_d = pseudo_dens;
375  auto phis_d = phis;
376 
377  // Input/Output data structures
378  //=======================================================
379  auto horiz_wind_d = horiz_wind;
380  auto cldfrac_liq_d = cldfrac_liq;
381  auto tke_d = tke;
382  auto qc_d = qc;
383 
384  // Enforce the correct grid heights and density
385  //=======================================================
386  auto dz_d = m_buffer.dz;
387 
388  // Subsidence pointer to device vector data
389  Real* w_sub = m_w_subsid;
390 
391  int nlay = m_num_layers;
392  Real dz = m_geom.CellSize(2);
393  bool moist = (m_cons->nComp() > RhoQ1_comp);
394  auto ProbLoArr = m_geom.ProbLoArray();
395 
396  auto domain = m_geom.Domain();
397  int ilo = domain.smallEnd(0);
398  int ihi = domain.bigEnd(0);
399  int jlo = domain.smallEnd(1);
400  int jhi = domain.bigEnd(1);
401 
402  for (MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
403  // NOTE: Grown box to get ghost cells in views
404  const auto& gbx = mfi.tilebox(IntVect(0,0,0),IntVect(1,1,0));
405  const int nx = gbx.length(0);
406  const int imin = gbx.smallEnd(0);
407  const int jmin = gbx.smallEnd(1);
408  const int kmax = gbx.bigEnd(2);
409  const int offset = m_col_offsets[mfi.index()];
410 
411  const Array4<const Real>& cons_arr = m_cons->const_array(mfi);
412 
413  const Array4<const Real>& u_arr = m_xvel->const_array(mfi);
414  const Array4<const Real>& v_arr = m_yvel->const_array(mfi);
415  const Array4<const Real>& w_arr = m_zvel->const_array(mfi);
416 
417  const Array4<const Real>& t13_arr = m_tau13->const_array(mfi);
418  const Array4<const Real>& t23_arr = m_tau23->const_array(mfi);
419  const Array4<const Real>& hfx3_arr = m_hfx3->const_array(mfi);
420  const Array4<const Real>& qfx3_arr = m_qfx3->const_array(mfi);
421 
422  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
423  Array4<const Real>{};
424  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
425  {
426  // NOTE: k gets permuted with ilay
427  // map [i,j,k] 0-based to [icol, ilay] 0-based
428  const int icol = (j-jmin)*nx + (i-imin) + offset;
429  const int ilay = kmax - k;
430  const int ilayi = kmax + 1 - k;
431 
432  // EOS input (at CC)
433  Real r = cons_arr(i,j,k,Rho_comp);
434  Real rt = cons_arr(i,j,k,RhoTheta_comp);
435  Real qv = (moist) ? cons_arr(i,j,k,RhoQ1_comp)/r : Real(0.);
436  Real qc = (moist) ? cons_arr(i,j,k,RhoQ2_comp)/r : Real(0.);
437 
438  // EOS avg to z-face
439  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
440  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
441  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : Real(0.);
442  Real rt_avg = Real(0.5) * (rt + rt_lo);
443  Real qv_avg = Real(0.5) * (qv + qv_lo);
444 
445  // Delta z
446  Real delz = (z_arr) ? Real(0.25) * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
447  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
448  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
449  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
450 
451  // W at cc (cannot be 0?; inspection of shoc code...)
452  Real w_cc = Real(0.5) * (w_arr(i,j,k) + w_arr(i,j,k+1));
453  w_cc += (w_sub) ? w_sub[k] : Real(0.);
454  Real w_limited = std::copysign(std::max(std::fabs(w_cc),Real(1.0e-6)),w_cc);
455 
456  // Input/Output data structures
457  //=======================================================
458  horiz_wind_d(icol,0,ilay) = Real(0.5) * (u_arr(i,j,k) + u_arr(i+1,j ,k));
459  horiz_wind_d(icol,1,ilay) = Real(0.5) * (v_arr(i,j,k) + v_arr(i ,j+1,k));
460  cldfrac_liq_d(icol,ilay) = (qc>Real(0.)) ? Real(1.) : Real(0.);
461  tke_d(icol,ilay) = std::max(cons_arr(i,j,k,RhoKE_comp)/r, Real(0.));
462  qc_d(icol,ilay) = qc;
463 
464  // Interface data structures
465  //=======================================================
466  // eamxx_common_physics_functions_impl.hpp: calculate_vertical_velocity
467  omega_d(icol,ilay) = -w_limited * r * CONST_GRAV;
468  if (k==0) {
469  int ii = std::min(std::max(i,ilo),ihi);
470  int jj = std::min(std::max(j,jlo),jhi);
471 
472  surf_mom_flux_d(icol,0) = Real(0.5) * (t13_arr(ii,jj,k) + t13_arr(ii+1,jj ,k));
473  surf_mom_flux_d(icol,1) = Real(0.5) * (t23_arr(ii,jj,k) + t23_arr(ii ,jj+1,k));
474  // No unit conversion to W/m^2 (ERF_ShocInterface.H L224)
475  surf_sens_flux_d(icol) = hfx3_arr(ii,jj,k);
476  surf_evap_d(icol) = (moist) ? qfx3_arr(ii,jj,k) : Real(0.);
477  // Back out the drag coeff
478  Real wsp = sqrt( horiz_wind_d(icol,0,ilay)[0]*horiz_wind_d(icol,0,ilay)[0]
479  + horiz_wind_d(icol,1,ilay)[0]*horiz_wind_d(icol,1,ilay)[0] );
480  surf_drag_coeff_tms_d(icol) = surf_mom_flux_d(icol,0) /
481  (-r * wsp * horiz_wind_d(icol,0,ilay)[0]);
482  }
483  T_mid_d(icol,ilay) = getTgivenRandRTh(r, rt, qv);
484  qv_d(icol,ilay) = qv;
485 
486  // Input data structures
487  //=======================================================
488  p_mid_d(icol,ilay) = getPgivenRTh(rt, qv);
489  p_int_d(icol,ilayi) = getPgivenRTh(rt_avg, qv_avg);
490  // eamxx_common_physics_functions_impl.hpp: calculate_density
491  pseudo_dens_d(icol,ilay) = r * CONST_GRAV * delz;
492  // Enforce the grid spacing
493  dz_d(icol,ilay) = delz;
494  // Surface geopotential
495  if (k==0) {
496  Real z = (z_arr) ? Real(0.125) * ( (z_arr(i ,j ,k+1) + z_arr(i ,j ,k))
497  + (z_arr(i+1,j ,k+1) + z_arr(i+1,j ,k))
498  + (z_arr(i ,j+1,k+1) + z_arr(i ,j+1,k))
499  + (z_arr(i+1,j+1,k+1) + z_arr(i+1,j+1,k)) ) : ProbLoArr[2];
500  phis_d(icol) = CONST_GRAV * z;
501  }
502 
503  if (ilay==0) {
504  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
505  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
506  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,Real(0.)) : Real(0.);
507  rt_avg = Real(0.5) * (rt + rt_hi);
508  qv_avg = Real(0.5) * (qv + qv_hi);
509  p_int_d(icol,0) = getPgivenRTh(rt_avg, qv_avg);
510  }
511  });
512  }
513 }
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:31
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
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
amrex::MultiFab * m_z_phys
Definition: ERF_ShocInterface.H:611
amrex::MultiFab * m_qfx3
Definition: ERF_ShocInterface.H:605
amrex::MultiFab * m_zvel
Definition: ERF_ShocInterface.H:596
amrex::MultiFab * m_tau23
Definition: ERF_ShocInterface.H:603
amrex::MultiFab * m_tau13
Definition: ERF_ShocInterface.H:602
amrex::MultiFab * m_hfx3
Definition: ERF_ShocInterface.H:604
amrex::Real * m_w_subsid
Definition: ERF_ShocInterface.H:599
Here is the call graph for this function:

◆ name()

std::string SHOCInterface::name ( ) const
inline
117 { return "shoc"; }

◆ requested_buffer_size_in_bytes()

size_t SHOCInterface::requested_buffer_size_in_bytes ( ) const
protected
761 {
762  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
763 
764  const int nlev_packs = ekat::npack<Spack>(m_num_layers);
765  const int nlevi_packs = ekat::npack<Spack>(m_num_layers+1);
766  const int num_tracer_packs = ekat::npack<Spack>(m_num_tracers);
767 
768  // Number of Reals needed by local views in the interface
769  const size_t interface_request = Buffer::num_1d_scalar_ncol*m_num_cols*sizeof(Real) +
770  Buffer::num_1d_scalar_nlev*nlev_packs*sizeof(Spack) +
771  Buffer::num_2d_vector_mid*m_num_cols*nlev_packs*sizeof(Spack) +
772  Buffer::num_2d_vector_int*m_num_cols*nlevi_packs*sizeof(Spack) +
773  Buffer::num_2d_vector_tr*m_num_cols*num_tracer_packs*sizeof(Spack);
774 
775  // Number of Reals needed by the WorkspaceManager passed to shoc_main
776  const auto policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
777  const int n_wind_slots = ekat::npack<Spack>(m_num_vel_comp)*Spack::n;
778  const int n_trac_slots = ekat::npack<Spack>(m_num_tracers) *Spack::n;
779  const size_t wsm_request = WSM::get_total_bytes_needed(nlevi_packs, 14+(n_wind_slots+n_trac_slots), policy);
780 
781  return ( (interface_request + wsm_request)/sizeof(Real) );
782 }
static constexpr int num_2d_vector_tr
Definition: ERF_ShocInterface.H:472
static constexpr int num_1d_scalar_nlev
Definition: ERF_ShocInterface.H:464

◆ run_impl()

void SHOCInterface::run_impl ( const Real  dt)
1022 {
1023  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
1024 
1025  EKAT_REQUIRE_MSG (dt<=300,
1026  "Error! SHOC is intended to run with a timestep no longer than 5 minutes.\n"
1027  " Please, reduce timestep (perhaps increasing subcycling iterations).\n");
1028 
1029  const auto nlev_packs = ekat::npack<Spack>(m_num_layers);
1030  const auto scan_policy = TPF::get_thread_range_parallel_scan_team_policy(m_num_cols, nlev_packs);
1031  const auto default_policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
1032 
1033  // Preprocessing of SHOC inputs. Kernel contains a parallel_scan,
1034  // so a special TeamPolicy is required.
1035  Kokkos::parallel_for("shoc_preprocess",
1036  scan_policy,
1037  shoc_preprocess);
1038  Kokkos::fence();
1039 
1040  auto wtracer_sfc = shoc_preprocess.wtracer_sfc;
1041  Kokkos::deep_copy(wtracer_sfc, 0);
1042 
1043  if (apply_tms) {
1045  }
1046 
1047  if (check_flux_state) {
1049  }
1050 
1051  // For now set the host timestep to the shoc timestep. This forces
1052  // number of SHOC timesteps (nadv) to be one
1053  // TODO: input parameter?
1054  hdtime = dt;
1055  m_nadv = std::max(static_cast<int>(round(hdtime/dt)),1);
1056 
1057  // Reset internal WSM variables.
1058  workspace_mgr.reset_internals();
1059 
1060  // Run shoc main
1061  SHF::shoc_main(m_num_cols, m_num_layers, m_num_layers+1, m_npbl, m_nadv, m_num_tracers, dt,
1063 #ifdef SCREAM_SHOC_SMALL_KERNELS
1064  , temporaries
1065 #endif
1066  );
1067 
1068  // Postprocessing of SHOC outputs
1069  Kokkos::parallel_for("shoc_postprocess",
1070  default_policy,
1072  Kokkos::fence();
1073 
1074  // Extra SHOC output diagnostics
1075  if (runtime_options.extra_diags) {
1076  Kokkos::deep_copy(shoc_mix,history_output.shoc_mix);
1077  Kokkos::deep_copy(brunt,history_output.brunt);
1078  Kokkos::deep_copy(w3,history_output.w3);
1079  Kokkos::deep_copy(isotropy,history_output.isotropy);
1080  Kokkos::deep_copy(wthl_sec,history_output.wthl_sec);
1081  Kokkos::deep_copy(wqw_sec,history_output.wqw_sec);
1082  Kokkos::deep_copy(uw_sec,history_output.uw_sec);
1083  Kokkos::deep_copy(vw_sec,history_output.vw_sec);
1084  Kokkos::deep_copy(qw_sec,history_output.qw_sec);
1085  Kokkos::deep_copy(thl_sec,history_output.thl_sec);
1086  }
1087 }
Int m_nadv
Definition: ERF_ShocInterface.H:573
void check_flux_state_consistency(const double dt)
Definition: ERF_ShocInterface.cpp:1126
Int hdtime
Definition: ERF_ShocInterface.H:576
void apply_turbulent_mountain_stress()
Definition: ERF_ShocInterface.cpp:1104

◆ set_computed_group_impl()

void SHOCInterface::set_computed_group_impl ( )
protected

◆ set_diff_stresses()

void SHOCInterface::set_diff_stresses ( )
672 {
673  for (MFIter mfi(*m_hfx3); mfi.isValid(); ++mfi) {
674  const auto& vbx_cc = mfi.validbox();
675  const auto& vbx_xz = convert(vbx_cc,IntVect(1,0,1));
676  const auto& vbx_yz = convert(vbx_cc,IntVect(0,1,1));
677 
678  const Array4<Real>& hfx_arr = m_hfx3->array(mfi);
679  const Array4<Real>& qfx_arr = m_qfx3->array(mfi);
680 
681  const Array4<Real>& t13_arr = m_tau13->array(mfi);
682  const Array4<Real>& t23_arr = m_tau23->array(mfi);
683 
684  ParallelFor(vbx_cc, vbx_xz, vbx_yz,
685  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
686  {
687  hfx_arr(i,j,k) = Real(0.);
688  qfx_arr(i,j,k) = Real(0.);
689  },
690  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
691  {
692  t13_arr(i,j,k) = Real(0.);
693  },
694  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
695  {
696  t23_arr(i,j,k) = Real(0.);
697  });
698  }
699 }
Here is the call graph for this function:

◆ set_eddy_diffs()

void SHOCInterface::set_eddy_diffs ( )
612 {
613  //
614  // Expose for device capture
615  //
616 
617  // Input/Output data structures
618  //=======================================================
619  auto tk_d = tk;
620 
621  // NOTE: Loop on grown box to fill ghost cells but limit
622  // to valid box where views are defined.
623  for (MFIter mfi(*m_mu); mfi.isValid(); ++mfi) {
624  const auto& gbx_cc = mfi.growntilebox();
625  const auto& vbx_cc = mfi.validbox();
626 
627  // NOTE: Grown box only for mapping
628  const auto& gbx = mfi.tilebox(IntVect(0,0,0),IntVect(1,1,0));
629  const int nx = gbx.length(0);
630  const int imin = gbx.smallEnd(0);
631  const int jmin = gbx.smallEnd(1);
632  const int kmax = gbx.bigEnd(2);
633  const int offset = m_col_offsets[mfi.index()];
634 
635  // Limiting to validbox
636  const int iminv = vbx_cc.smallEnd(0);
637  const int imaxv = vbx_cc.bigEnd(0);
638  const int jminv = vbx_cc.smallEnd(1);
639  const int jmaxv = vbx_cc.bigEnd(1);
640  const int kminv = vbx_cc.smallEnd(2);
641  const int kmaxv = vbx_cc.bigEnd(2);
642 
643  const Array4<Real>& mu_arr = m_mu->array(mfi);
644 
645  ParallelFor(gbx_cc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
646  {
647  // Limiting
648  int ii = std::min(std::max(i,iminv),imaxv);
649  int jj = std::min(std::max(j,jminv),jmaxv);
650  int kk = std::min(std::max(k,kminv),kmaxv);
651 
652  // NOTE: k gets permuted with ilay
653  // map [i,j,k] 0-based to [icol, ilay] 0-based
654  const int icol = (jj-jmin)*nx + (ii-imin) + offset;
655  const int ilay = kmax - kk;
656 
657  // NOTE: Set mom_v for tau_33, all other vertical comps are 0
658  mu_arr(i,j,k,EddyDiff::Mom_v) = tk_d(icol,ilay)[0];
659  mu_arr(i,j,k,EddyDiff::Theta_v) = Real(0.);
660  mu_arr(i,j,k,EddyDiff::KE_v) = Real(0.);
661  mu_arr(i,j,k,EddyDiff::Q_v) = Real(0.);
662  });
663  }
664 
665  // Correct the internal ghost cells that have foextrap
666  m_mu->FillBoundary(m_geom.periodicity());
667 }
amrex::MultiFab * m_mu
Definition: ERF_ShocInterface.H:608
@ Theta_v
Definition: ERF_IndexDefines.H:194
@ Q_v
Definition: ERF_IndexDefines.H:197
@ Mom_v
Definition: ERF_IndexDefines.H:193
@ KE_v
Definition: ERF_IndexDefines.H:195
Here is the call graph for this function:

◆ set_grids()

void SHOCInterface::set_grids ( int &  level,
const amrex::BoxArray &  ba,
amrex::Geometry &  geom,
amrex::MultiFab *  cons,
amrex::MultiFab *  xvel,
amrex::MultiFab *  yvel,
amrex::MultiFab *  zvel,
amrex::Real w_subsid,
amrex::MultiFab *  tau13,
amrex::MultiFab *  tau23,
amrex::MultiFab *  hfx3,
amrex::MultiFab *  qfx3,
amrex::MultiFab *  eddyDiffs,
amrex::MultiFab *  z_phys 
)
135 {
136  // Set data members that may change
137  m_lev = level;
138  m_geom = geom;
139  m_cons = cons;
140  m_xvel = xvel;
141  m_yvel = yvel;
142  m_zvel = zvel;
143  m_w_subsid = w_subsid;
144  m_tau13 = tau13;
145  m_tau23 = tau23;
146  m_hfx3 = hfx3;
147  m_qfx3 = qfx3;
148  m_mu = eddyDiffs;
149  m_z_phys = z_phys;
150 
151  // Ensure the boxes span klo -> khi
152  int klo = geom.Domain().smallEnd(2);
153  int khi = geom.Domain().bigEnd(2);
154 
155  // Reset vector of offsets for columnar data
156  m_num_layers = geom.Domain().length(2);
157 
158  int num_cols = 0;
159  m_col_offsets.clear();
160  m_col_offsets.resize(int(ba.size()));
161  for (MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
162  // NOTE: Get lateral ghost cells for CC <--> FC
163  const auto& gbx = mfi.tilebox(IntVect(0,0,0),IntVect(1,1,0));
164  AMREX_ALWAYS_ASSERT_WITH_MESSAGE((klo == gbx.smallEnd(2)) &&
165  (khi == gbx.bigEnd(2)),
166  "Vertical decomposition with shoc is not allowed.");
167  int nx = gbx.length(0);
168  int ny = gbx.length(1);
169  m_col_offsets[mfi.index()] = num_cols;
170  num_cols += nx * ny;
171  }
172 
173  // Resize the Kokkos variables that persist in memory
174  if (num_cols != m_num_cols) {
176  tk = view_2d();
177  sgs_buoy_flux = view_2d("sgs_buoy_flux", num_cols, m_num_layers);
178  tk = view_2d("eddy_diff_mom", num_cols, m_num_layers);
179  }
180  m_num_cols = num_cols;
181 
182  // Allocate the tendency MultiFabs
183  c_tend.define(m_cons->boxArray(), m_cons->DistributionMap(), m_cons->nComp(), 0);
184  u_tend.define(m_xvel->boxArray(), m_xvel->DistributionMap(), m_xvel->nComp(), 0);
185  v_tend.define(m_yvel->boxArray(), m_yvel->DistributionMap(), m_yvel->nComp(), 0);
186 
187  // Allocate the buffer arrays in ERF
188  alloc_buffers();
189 
190  // Allocate the m_buffer struct
191  init_buffers();
192 
193  // Fill the KOKKOS Views from AMReX MFs
195 }
@ tau23
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
void alloc_buffers()
Definition: ERF_ShocInterface.cpp:199
void init_buffers()
Definition: ERF_ShocInterface.cpp:786
int m_lev
Definition: ERF_ShocInterface.H:579
void mf_to_kokkos_buffers()
Definition: ERF_ShocInterface.cpp:348
@ xvel
Definition: ERF_IndexDefines.H:159
@ cons
Definition: ERF_IndexDefines.H:158
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160

Member Data Documentation

◆ apply_tms

bool SHOCInterface::apply_tms = false
protected

◆ brunt

view_2d SHOCInterface::brunt
protected

◆ c_tend

amrex::MultiFab SHOCInterface::c_tend
protected

◆ check_flux_state

bool SHOCInterface::check_flux_state = false
protected

◆ cldfrac_liq

view_2d SHOCInterface::cldfrac_liq
protected

◆ cldfrac_liq_prev

view_2d SHOCInterface::cldfrac_liq_prev
protected

◆ column_conservation_check

bool SHOCInterface::column_conservation_check = false
protected

◆ extra_shoc_diags

bool SHOCInterface::extra_shoc_diags = false
protected

◆ hdtime

Int SHOCInterface::hdtime
protected

◆ heat_flux

view_1d SHOCInterface::heat_flux
protected

◆ history_output

SHF::SHOCHistoryOutput SHOCInterface::history_output
protected

◆ horiz_wind

view_3d SHOCInterface::horiz_wind
protected

◆ ice_flux

view_1d SHOCInterface::ice_flux
protected

◆ input

SHF::SHOCInput SHOCInterface::input
protected

◆ input_output

SHF::SHOCInputOutput SHOCInterface::input_output
protected

◆ inv_qc_relvar

view_2d SHOCInterface::inv_qc_relvar
protected

◆ isotropy

view_2d SHOCInterface::isotropy
protected

◆ m_ba

amrex::BoxArray SHOCInterface::m_ba
protected

◆ m_buffer

Buffer SHOCInterface::m_buffer
protected

◆ m_col_offsets

amrex::Vector<int> SHOCInterface::m_col_offsets
protected

◆ m_cons

amrex::MultiFab* SHOCInterface::m_cons = nullptr
protected

◆ m_first_step

bool SHOCInterface::m_first_step = true
protected

◆ m_geom

amrex::Geometry SHOCInterface::m_geom
protected

◆ m_hfx3

amrex::MultiFab* SHOCInterface::m_hfx3 = nullptr
protected

◆ m_lev

int SHOCInterface::m_lev
protected

◆ m_mu

amrex::MultiFab* SHOCInterface::m_mu = nullptr
protected

◆ m_nadv

Int SHOCInterface::m_nadv
protected

◆ m_npbl

Int SHOCInterface::m_npbl
protected

◆ m_num_cols

Int SHOCInterface::m_num_cols = 0
protected

◆ m_num_layers

Int SHOCInterface::m_num_layers = 0
protected

◆ m_num_tracers

Int SHOCInterface::m_num_tracers = 3
protected

◆ m_num_vel_comp

Int SHOCInterface::m_num_vel_comp = 2
protected

◆ m_qfx3

amrex::MultiFab* SHOCInterface::m_qfx3 = nullptr
protected

◆ m_step

int SHOCInterface::m_step
protected

◆ m_tau13

amrex::MultiFab* SHOCInterface::m_tau13 = nullptr
protected

◆ m_tau23

amrex::MultiFab* SHOCInterface::m_tau23 = nullptr
protected

◆ m_w_subsid

amrex::Real* SHOCInterface::m_w_subsid = nullptr
protected

◆ m_xvel

amrex::MultiFab* SHOCInterface::m_xvel = nullptr
protected

◆ m_yvel

amrex::MultiFab* SHOCInterface::m_yvel = nullptr
protected

◆ m_z_phys

amrex::MultiFab* SHOCInterface::m_z_phys = nullptr
protected

◆ m_zvel

amrex::MultiFab* SHOCInterface::m_zvel = nullptr
protected

◆ obklen

view_1d SHOCInterface::obklen
protected

◆ omega

view_2d SHOCInterface::omega
protected

◆ output

SHF::SHOCOutput SHOCInterface::output
protected

◆ p_int

view_2d SHOCInterface::p_int
protected

◆ p_mid

view_2d SHOCInterface::p_mid
protected

◆ pblh

view_1d SHOCInterface::pblh
protected

◆ phis

view_1d SHOCInterface::phis
protected

◆ pseudo_dens

view_2d SHOCInterface::pseudo_dens
protected

◆ qc

view_2d SHOCInterface::qc
protected

◆ qtracers

view_3d_strided SHOCInterface::qtracers
protected

◆ qv

view_2d SHOCInterface::qv
protected

◆ qw_sec

view_2d SHOCInterface::qw_sec
protected

◆ runtime_options

SHF::SHOCRuntime SHOCInterface::runtime_options
protected

◆ sgs_buoy_flux

view_2d SHOCInterface::sgs_buoy_flux
protected

◆ shoc_cond

view_2d SHOCInterface::shoc_cond
protected

◆ shoc_evap

view_2d SHOCInterface::shoc_evap
protected

◆ shoc_mix

view_2d SHOCInterface::shoc_mix
protected

◆ shoc_postprocess

SHOCPostprocess SHOCInterface::shoc_postprocess
protected

◆ shoc_preprocess

SHOCPreprocess SHOCInterface::shoc_preprocess
protected

◆ surf_drag_coeff_tms

view_1d SHOCInterface::surf_drag_coeff_tms
protected

◆ surf_evap

view_1d SHOCInterface::surf_evap
protected

◆ surf_mom_flux

sview_2d SHOCInterface::surf_mom_flux
protected

◆ surf_sens_flux

view_1d SHOCInterface::surf_sens_flux
protected

◆ T_mid

view_2d SHOCInterface::T_mid
protected

◆ thl_sec

view_2d SHOCInterface::thl_sec
protected

◆ tk

view_2d SHOCInterface::tk
protected

◆ tke

view_2d SHOCInterface::tke
protected

◆ tkh

view_2d SHOCInterface::tkh
protected

◆ tot_buff_view

view_1d SHOCInterface::tot_buff_view
protected

◆ u_tend

amrex::MultiFab SHOCInterface::u_tend
protected

◆ ustar

view_1d SHOCInterface::ustar
protected

◆ uw_sec

view_2d SHOCInterface::uw_sec
protected

◆ v_tend

amrex::MultiFab SHOCInterface::v_tend
protected

◆ vapor_flux

view_1d SHOCInterface::vapor_flux
protected

◆ vw_sec

view_2d SHOCInterface::vw_sec
protected

◆ w3

view_2d SHOCInterface::w3
protected

◆ w_sec

view_2d SHOCInterface::w_sec
protected

◆ water_flux

view_1d SHOCInterface::water_flux
protected

◆ workspace_mgr

ekat::WorkspaceManager<Spack, KT::Device> SHOCInterface::workspace_mgr
protected

◆ wqw_sec

view_2d SHOCInterface::wqw_sec
protected

◆ wthl_sec

view_2d SHOCInterface::wthl_sec
protected

The documentation for this class was generated from the following files: