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 *z_phys)
 
void initialize_impl ()
 
void run_impl (const Real dt)
 
void finalize_impl ()
 
void alloc_buffers ()
 
void dealloc_buffers ()
 
void mf_to_kokkos_buffers ()
 
void kokkos_buffers_to_mf ()
 
std::string name () const
 

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
 
Int m_num_layers
 
Int m_npbl
 
Int m_nadv
 
Int m_num_tracers = 1
 
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_z_phys = nullptr
 
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
 
bool apply_tms = false
 
bool check_flux_state = false
 
view_1d phis
 
view_1d wpthlp_sfc
 
view_1d wprtp_sfc
 
view_1d upwp_sfc
 
view_1d vpwp_sfc
 
view_1d surf_sens_flux
 
view_1d surf_evap
 
sview_2d surf_mom_flux
 
view_2d T_mid
 
view_2d p_mid
 
view_2d p_int
 
view_2d dens
 
view_2d omega
 
view_2d qv
 
view_2d qc
 
view_2d qc_copy
 
view_2d z_mid
 
view_2d z_int
 
view_2d shoc_s
 
view_2d tke
 
view_2d tke_copy
 
view_2d rrho
 
view_2d rrho_i
 
view_2d thv
 
view_2d dz
 
view_2d dse
 
view_2d zt_grid
 
view_2d zi_grid
 
view_2d wtracer_sfc
 
view_2d wm_zt
 
view_2d inv_exner
 
view_2d thlm
 
view_2d qw
 
view_2d cloud_frac
 
view_2d cldfrac_liq
 
view_2d cldfrac_liq_prev
 
view_3d_strided qtracers
 
view_2d shoc_ql2
 
view_2d inv_qc_relvar
 
view_2d sgs_buoy_flux
 
view_2d tk
 
view_3d horiz_wind
 
view_1d pblh
 
view_1d ustar
 
view_1d obklen
 
view_2d tkh
 
view_2d shoc_mix
 
view_2d w_sec
 
view_2d thl_sec
 
view_2d qw_sec
 
view_2d qwthl_sec
 
view_2d wthl_sec
 
view_2d wqw_sec
 
view_2d wtke_sec
 
view_2d uw_sec
 
view_2d vw_sec
 
view_2d w3
 
view_2d wqls_sec
 
view_2d brunt
 
view_2d isotropy
 
view_2d shoc_cond
 
view_2d shoc_evap
 
view_1d se_b
 
view_1d ke_b
 
view_1d wv_b
 
view_1d wl_b
 
view_1d se_a
 
view_1d ke_a
 
view_1d wv_a
 
view_1d wl_a
 
view_1d kbfs
 
view_1d ustar2
 
view_1d wstar
 
view_1d surf_drag_coeff_tms
 
view_2d rho_zt
 
view_2d shoc_qv
 
view_2d tabs
 
view_2d dz_zt
 
view_2d dz_zi
 

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 ekat::KokkosTypes< KokkosDefaultDevice >::template view_2d< Real >
 
using sview_2d_const = typename ekat::KokkosTypes< KokkosDefaultDevice >::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 ekat::KokkosTypes<KokkosDefaultDevice>::template view_2d<Real>
private

◆ sview_2d_const

using SHOCInterface::sview_2d_const = typename ekat::KokkosTypes<KokkosDefaultDevice>::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  bool def_check_flux_state = false;
13  bool def_enable_column_conservation = false;
14  // Extra SHOC diagnostics
15  bool def_extra_diags = false;
16  // Turn off SGS variability in SHOC, effectively reducing it to a 1.5 TKE closure?
17  bool def_shoc_1p5tke = false;
18  // Minimum value of stability correction
19  Real def_lambda_low = 0.001;
20  // Maximum value of stability correction
21  Real def_lambda_high = 0.08;
22  // Slope of change from lambda_low to lambda_high
23  Real def_lambda_slope = 0.08;
24  // stability threshold for which to apply more stability correction
25  Real def_lambda_thresh = 0.02;
26  // Temperature variance tuning factor
27  Real def_thl2tune = 1.0;
28  // Moisture variance tuning factor
29  Real def_qw2tune = 1.0;
30  // Temperature moisture covariance
31  Real def_qwthl2tune = 1.0;
32  // Vertical velocity variance
33  Real def_w2tune = 1.0;
34  // Length scale factor
35  Real def_length_fac = 0.5;
36  // Third moment vertical velocity damping factor
37  Real def_c_diag_3rd_mom = 7.0;
38  // Eddy diffusivity coefficient for heat
39  Real def_coeff_kh = 0.1;
40  // Eddy diffusivity coefficient for momentum
41  Real def_coeff_km = 0.1;
42 
43  // From E3SM/components/eamxx/src/control/atmospheric_driver.cpp
44  bool def_apply_tms = true;
45 
46  runtime_options.lambda_low = def_lambda_low;
47  runtime_options.lambda_high = def_lambda_high;
48  runtime_options.lambda_slope = def_lambda_slope;
49  runtime_options.lambda_thresh = def_lambda_thresh;
50 
51  runtime_options.thl2tune = def_thl2tune;
52  runtime_options.qwthl2tune = def_qwthl2tune;
53  runtime_options.qw2tune = def_qw2tune;
54  runtime_options.w2tune = def_w2tune;
55 
56  runtime_options.length_fac = def_length_fac;
57  runtime_options.c_diag_3rd_mom = def_c_diag_3rd_mom;
58  runtime_options.Ckh = def_coeff_kh;
59  runtime_options.Ckm = def_coeff_km;
60  runtime_options.shoc_1p5tke = def_shoc_1p5tke;
61  runtime_options.extra_diags = def_extra_diags;
62 
63  // Construct parser object for following reads
64  ParmParse pp("erf.shoc");
65 
66  // Parse runtime inputs at start up
67  pp.query("lambda_low" , runtime_options.lambda_low );
68  pp.query("lambda_high" , runtime_options.lambda_high );
69  pp.query("lambda_slope" , runtime_options.lambda_slope );
70  pp.query("lambda_thresh" , runtime_options.lambda_thresh );
71  pp.query("thl2tune" , runtime_options.thl2tune );
72  pp.query("qw2tune" , runtime_options.qw2tune );
73  pp.query("qwthl2tune" , runtime_options.qwthl2tune );
74  pp.query("w2tune" , runtime_options.w2tune );
75  pp.query("length_fac" , runtime_options.length_fac );
76  pp.query("c_diag_3rd_mom" , runtime_options.c_diag_3rd_mom);
77  pp.query("coeff_kh" , runtime_options.Ckh );
78  pp.query("coeff_km" , runtime_options.Ckm );
79  pp.query("shoc_1p5tke" , runtime_options.shoc_1p5tke );
80  pp.query("extra_shoc_diags", runtime_options.extra_diags );
81 
82  // Set to default but allow us to change it through the inputs file
83  apply_tms = def_apply_tms;
84  pp.query("apply_tms", apply_tms);
85 
86  pp.query("check_flux_state", check_flux_state);
87 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:230
amrex::Real Real
Definition: ERF_ShocInterface.H:19
bool check_flux_state
Definition: ERF_ShocInterface.H:588
bool apply_tms
Definition: ERF_ShocInterface.H:587
SHF::SHOCRuntime runtime_options
Definition: ERF_ShocInterface.H:573
Here is the call graph for this function:

Member Function Documentation

◆ alloc_buffers()

void SHOCInterface::alloc_buffers ( )
155 {
156  // SHOCPreprocess data structures
157  //=======================================================
158  phis = view_1d("Phis" , m_num_cols);
159  wpthlp_sfc = view_1d("Wpthlp_sfc" , m_num_cols);
160  wprtp_sfc = view_1d("Wprtp_sfc" , m_num_cols);
161  upwp_sfc = view_1d("Upwp_sfc" , m_num_cols);
162  vpwp_sfc = view_1d("Vpwp_sfc" , m_num_cols);
163  surf_sens_flux = view_1d("Sfc sens flux" , m_num_cols);
164  surf_evap = view_1d("Sfc latent flux", m_num_cols);
165  surf_mom_flux = sview_2d("Sfc mom flux" , m_num_cols, m_num_vel_comp);
166  wtracer_sfc = view_2d("Wtracer_sfc" , m_num_cols, m_num_tracers );
167 
168  T_mid = view_2d("T_mid" , m_num_cols, m_num_layers );
169  p_mid = view_2d("P_mid" , m_num_cols, m_num_layers );
170  p_int = view_2d("P_int" , m_num_cols, m_num_layers+1);
171  dens = view_2d("Density" , m_num_cols, m_num_layers );
172  omega = view_2d("Omega" , m_num_cols, m_num_layers );
173  qv = view_2d("Qv" , m_num_cols, m_num_layers );
174  qc = view_2d("Qc" , m_num_cols, m_num_layers );
175  qc_copy = view_2d("Qc_copy" , m_num_cols, m_num_layers );
176  z_mid = view_2d("Z_mid" , m_num_cols, m_num_layers );
177  z_int = view_2d("Z_int" , m_num_cols, m_num_layers+1);
178  shoc_s = view_2d("Shoc_s" , m_num_cols, m_num_layers );
179  tke = view_2d("Tke" , m_num_cols, m_num_layers );
180  tke_copy = view_2d("Tke_copy" , m_num_cols, m_num_layers );
181  rrho = view_2d("Rrho" , m_num_cols, m_num_layers );
182  rrho_i = view_2d("Rrho_i" , m_num_cols, m_num_layers );
183  thv = view_2d("Thv" , m_num_cols, m_num_layers );
184  dz = view_2d("dz" , m_num_cols, m_num_layers );
185  dse = view_2d("dse" , m_num_cols, m_num_layers );
186  zt_grid = view_2d("Zt_grid" , m_num_cols, m_num_layers );
187  zi_grid = view_2d("Zi_grid" , m_num_cols, m_num_layers );
188  wm_zt = view_2d("wm_zt" , m_num_cols, m_num_layers );
189  inv_exner = view_2d("Inv_exner" , m_num_cols, m_num_layers );
190  thlm = view_2d("Thlm" , m_num_cols, m_num_layers );
191  qw = view_2d("Qw" , m_num_cols, m_num_layers );
192  cloud_frac = view_2d("Cld_frac" , m_num_cols, m_num_layers );
193  cldfrac_liq = view_2d("Cld_frac_liq" , m_num_cols, m_num_layers );
194  cldfrac_liq_prev = view_2d("cld_frac_liq_prev" , m_num_cols, m_num_layers );
195 
196  // NOTE: Use layoutright format
197  Kokkos::LayoutStride layout(m_num_cols , m_num_cols*m_num_layers, // stride for dim0
198  m_num_layers , m_num_layers, // stride for dim1
199  m_num_tracers, 1 // stride for dim2
200  );
201  qtracers = view_3d_strided("Qtracers" , layout);
202 
203 
204  // SHOCPostprocess data structures
205  //=======================================================
207  inv_qc_relvar = view_2d("inv_qc_relvar", m_num_cols, m_num_layers);
208 
209  // SHOC InputOutput data structures
210  //=======================================================
211  sgs_buoy_flux = view_2d("sgs_buoy_flux", m_num_cols, m_num_layers);
212  tk = view_2d("eddy_diff_mom", m_num_cols, m_num_layers);
213 
215 
216  // SHOC Output data structures
217  //=======================================================
218  pblh = view_1d("pbl_height", m_num_cols);
219  ustar = view_1d("ustar" , m_num_cols);
220  obklen = view_1d("obklen" , m_num_cols);
221 
222  tkh = view_2d("eddy_diff_heat", m_num_cols, m_num_layers);
223 
224  // SHOC HistoryOutput data structures
225  //=======================================================
226  shoc_mix = view_2d("shoc_mix" , m_num_cols, m_num_layers);
227  w_sec = view_2d("w_sec" , m_num_cols, m_num_layers);
228  thl_sec = view_2d("thl_sec" , m_num_cols, m_num_layers);
229  qw_sec = view_2d("qw_sec" , m_num_cols, m_num_layers);
230  qwthl_sec = view_2d("qwthl_sec", m_num_cols, m_num_layers);
231  wthl_sec = view_2d("wthl_sec" , m_num_cols, m_num_layers);
232  wqw_sec = view_2d("wqw_sec" , m_num_cols, m_num_layers);
233  wtke_sec = view_2d("wtke_sec" , m_num_cols, m_num_layers);
234  uw_sec = view_2d("uw_sec" , m_num_cols, m_num_layers);
235  vw_sec = view_2d("vw_sec" , m_num_cols, m_num_layers);
236  w3 = view_2d("w3" , m_num_cols, m_num_layers);
237  wqls_sec = view_2d("wqls_sec" , m_num_cols, m_num_layers);
238  brunt = view_2d("brunt" , m_num_cols, m_num_layers);
239  isotropy = view_2d("isotropy" , m_num_cols, m_num_layers);
240  shoc_cond = view_2d("shoc_cond", m_num_cols, m_num_layers);
241  shoc_evap = view_2d("shoc_evap", m_num_cols, m_num_layers);
242 
243  // SHOC Miscellaneous data structures
244  //=======================================================
245  // Small kernels
246  se_b = view_1d("se_b" , m_num_cols);
247  ke_b = view_1d("ke_b" , m_num_cols);
248  wv_b = view_1d("wv_b" , m_num_cols);
249  wl_b = view_1d("wl_b" , m_num_cols);
250  se_a = view_1d("se_a" , m_num_cols);
251  ke_a = view_1d("ke_a" , m_num_cols);
252  wv_a = view_1d("wv_a" , m_num_cols);
253  wl_a = view_1d("wl_a" , m_num_cols);
254  kbfs = view_1d("kbfs" , m_num_cols);
255  ustar2 = view_1d("ustar2", m_num_cols);
256  wstar = view_1d("wstar" , m_num_cols);
257  if (apply_tms) { surf_drag_coeff_tms = view_1d("surf_drag_coeff" , m_num_cols); }
258 
259  rho_zt = view_2d("rho_zt" , m_num_cols, m_num_layers);
260  shoc_qv = view_2d("shoc_qv" , m_num_cols, m_num_layers);
261  tabs = view_2d("tabs" , m_num_cols, m_num_layers);
262  dz_zt = view_2d("dz_zt" , m_num_cols, m_num_layers);
263  dz_zi = view_2d("dz_zi" , m_num_cols, m_num_layers+1);
264 }
view_2d tabs
Definition: ERF_ShocInterface.H:686
view_2d p_int
Definition: ERF_ShocInterface.H:603
view_2d wqw_sec
Definition: ERF_ShocInterface.H:658
view_1d surf_drag_coeff_tms
Definition: ERF_ShocInterface.H:683
view_2d dz_zi
Definition: ERF_ShocInterface.H:688
view_2d inv_exner
Definition: ERF_ShocInterface.H:623
view_1d wl_b
Definition: ERF_ShocInterface.H:675
view_1d wv_a
Definition: ERF_ShocInterface.H:678
view_2d wtke_sec
Definition: ERF_ShocInterface.H:659
view_2d thl_sec
Definition: ERF_ShocInterface.H:654
view_2d tkh
Definition: ERF_ShocInterface.H:648
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:627
Int m_num_layers
Definition: ERF_ShocInterface.H:537
view_2d vw_sec
Definition: ERF_ShocInterface.H:661
view_2d thv
Definition: ERF_ShocInterface.H:616
view_2d shoc_qv
Definition: ERF_ShocInterface.H:685
view_2d shoc_cond
Definition: ERF_ShocInterface.H:666
typename ekat::KokkosTypes< KokkosDefaultDevice >::template view_2d< Real > sview_2d
Definition: ERF_ShocInterface.H:51
view_1d vpwp_sfc
Definition: ERF_ShocInterface.H:596
view_2d brunt
Definition: ERF_ShocInterface.H:664
Int m_num_vel_comp
Definition: ERF_ShocInterface.H:541
view_2d dens
Definition: ERF_ShocInterface.H:604
view_1d surf_sens_flux
Definition: ERF_ShocInterface.H:597
view_1d pblh
Definition: ERF_ShocInterface.H:645
Int m_num_cols
Definition: ERF_ShocInterface.H:536
view_1d wstar
Definition: ERF_ShocInterface.H:682
view_1d se_b
Definition: ERF_ShocInterface.H:672
view_2d sgs_buoy_flux
Definition: ERF_ShocInterface.H:639
view_2d p_mid
Definition: ERF_ShocInterface.H:602
view_2d rho_zt
Definition: ERF_ShocInterface.H:684
view_2d uw_sec
Definition: ERF_ShocInterface.H:660
view_2d dse
Definition: ERF_ShocInterface.H:618
view_2d qc
Definition: ERF_ShocInterface.H:607
view_1d kbfs
Definition: ERF_ShocInterface.H:680
view_2d dz
Definition: ERF_ShocInterface.H:617
view_2d qwthl_sec
Definition: ERF_ShocInterface.H:656
view_1d ke_a
Definition: ERF_ShocInterface.H:677
view_1d wl_a
Definition: ERF_ShocInterface.H:679
view_2d shoc_mix
Definition: ERF_ShocInterface.H:652
view_2d qw_sec
Definition: ERF_ShocInterface.H:655
view_2d cldfrac_liq_prev
Definition: ERF_ShocInterface.H:628
view_2d qw
Definition: ERF_ShocInterface.H:625
view_2d cloud_frac
Definition: ERF_ShocInterface.H:626
view_2d omega
Definition: ERF_ShocInterface.H:605
typename SHF::view_3d_strided< Spack > view_3d_strided
Definition: ERF_ShocInterface.H:55
view_1d ustar2
Definition: ERF_ShocInterface.H:681
view_2d rrho_i
Definition: ERF_ShocInterface.H:615
view_1d obklen
Definition: ERF_ShocInterface.H:647
view_2d z_mid
Definition: ERF_ShocInterface.H:609
view_2d qc_copy
Definition: ERF_ShocInterface.H:608
view_2d rrho
Definition: ERF_ShocInterface.H:614
view_2d inv_qc_relvar
Definition: ERF_ShocInterface.H:635
view_1d surf_evap
Definition: ERF_ShocInterface.H:598
view_2d w_sec
Definition: ERF_ShocInterface.H:653
view_1d phis
Definition: ERF_ShocInterface.H:592
view_2d shoc_evap
Definition: ERF_ShocInterface.H:667
view_3d horiz_wind
Definition: ERF_ShocInterface.H:641
typename SHF::view_1d< Real > view_1d
Definition: ERF_ShocInterface.H:47
view_2d wqls_sec
Definition: ERF_ShocInterface.H:663
view_1d ke_b
Definition: ERF_ShocInterface.H:673
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:630
view_1d upwp_sfc
Definition: ERF_ShocInterface.H:595
view_1d ustar
Definition: ERF_ShocInterface.H:646
view_2d tke_copy
Definition: ERF_ShocInterface.H:613
view_2d wm_zt
Definition: ERF_ShocInterface.H:622
view_2d w3
Definition: ERF_ShocInterface.H:662
view_1d se_a
Definition: ERF_ShocInterface.H:676
view_1d wv_b
Definition: ERF_ShocInterface.H:674
Int m_num_tracers
Definition: ERF_ShocInterface.H:540
sview_2d surf_mom_flux
Definition: ERF_ShocInterface.H:599
typename SHF::view_2d< SHF::Spack > view_2d
Definition: ERF_ShocInterface.H:49
view_2d zi_grid
Definition: ERF_ShocInterface.H:620
view_2d thlm
Definition: ERF_ShocInterface.H:624
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:621
view_2d shoc_ql2
Definition: ERF_ShocInterface.H:634
view_2d tk
Definition: ERF_ShocInterface.H:640
view_2d qv
Definition: ERF_ShocInterface.H:606
view_2d wthl_sec
Definition: ERF_ShocInterface.H:657
view_2d tke
Definition: ERF_ShocInterface.H:612
view_2d zt_grid
Definition: ERF_ShocInterface.H:619
view_2d T_mid
Definition: ERF_ShocInterface.H:601
view_2d shoc_s
Definition: ERF_ShocInterface.H:611
view_2d z_int
Definition: ERF_ShocInterface.H:610
typename SHF::view_3d< Spack > view_3d
Definition: ERF_ShocInterface.H:53
view_2d isotropy
Definition: ERF_ShocInterface.H:665
view_2d dz_zt
Definition: ERF_ShocInterface.H:687
view_1d wpthlp_sfc
Definition: ERF_ShocInterface.H:593
view_1d wprtp_sfc
Definition: ERF_ShocInterface.H:594

◆ apply_turbulent_mountain_stress()

void SHOCInterface::apply_turbulent_mountain_stress ( )
protected
736 {
737  const int nlev_v = (m_num_layers-1)/Spack::n;
738  const int nlev_p = (m_num_layers-1)%Spack::n;
739  const int nlevi_v = m_num_layers/Spack::n;
740  const int nlevi_p = m_num_layers%Spack::n;
741 
742  Kokkos::parallel_for("apply_tms", KT::RangePolicy(0, m_num_cols), KOKKOS_LAMBDA (const int i)
743  {
744  upwp_sfc(i) -= surf_drag_coeff_tms(i)*horiz_wind(i,0,nlev_v)[nlev_p]/rrho_i(i,nlevi_v)[nlevi_p];
745  vpwp_sfc(i) -= surf_drag_coeff_tms(i)*horiz_wind(i,1,nlev_v)[nlev_p]/rrho_i(i,nlevi_v)[nlevi_p];
746  });
747 }

◆ check_flux_state_consistency()

void SHOCInterface::check_flux_state_consistency ( const double  dt)
protected
752 {
753  using PC = scream::physics::Constants<Real>;
754  using RU = ekat::ReductionUtils<KT::ExeSpace>;
755  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
756 
757  const Real gravit = PC::gravit;
758  const Real qmin = 1e-12; // minimum permitted constituent concentration (kg/kg)
759 
760  const auto& pseudo_density = dens;
761 
762  const auto nlevs = m_num_layers;
763  const auto nlev_packs = ekat::npack<Spack>(nlevs);
764  const auto last_pack_idx = (nlevs-1)/Spack::n;
765  const auto last_pack_entry = (nlevs-1)%Spack::n;
766  const auto policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
767  Kokkos::parallel_for("check_flux_state_consistency",
768  policy,
769  KOKKOS_LAMBDA (const KT::MemberType& team)
770  {
771  const auto i = team.league_rank();
772 
773  const auto& pseudo_density_i = ekat::subview(pseudo_density, i);
774  const auto& qv_i = ekat::subview(qv, i);
775 
776  // reciprocal of pseudo_density at the bottom layer
777  const auto rpdel = 1.0/pseudo_density_i(last_pack_idx)[last_pack_entry];
778 
779  // Check if the negative surface latent heat flux can exhaust
780  // the moisture in the lowest model level. If so, apply fixer.
781  const auto condition = surf_evap(i) - (qmin - qv_i(last_pack_idx)[last_pack_entry])/(dt*gravit*rpdel);
782  if (condition < 0) {
783  const auto cc = std::fabs(surf_evap(i)*dt*gravit);
784 
785  auto tracer_mass = [&](const int k)
786  {
787  return qv_i(k)*pseudo_density_i(k);
788  };
789  Real mm = RU::view_reduction(team, 0, nlevs, tracer_mass);
790 
791  EKAT_KERNEL_ASSERT_MSG(mm >= cc, "Error! Total mass of column vapor should be greater than mass of surf_evap.\n");
792 
793  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&](const int& k)
794  {
795  const auto adjust = cc*qv_i(k)*pseudo_density_i(k)/mm;
796  qv_i(k) = (qv_i(k)*pseudo_density_i(k) - adjust)/pseudo_density_i(k);
797  });
798 
799  surf_evap(i) = 0;
800  }
801  });
802 }
@ PC
Definition: ERF_IndexDefines.H:123

◆ dealloc_buffers()

void SHOCInterface::dealloc_buffers ( )
269 {
270  // SHOCPreprocess data structures
271  //=======================================================
272  phis = view_1d();
273  wpthlp_sfc = view_1d();
274  wprtp_sfc = view_1d();
275  upwp_sfc = view_1d();
276  vpwp_sfc = view_1d();
278  surf_evap = view_1d();
280  wtracer_sfc = view_2d();
281 
282  T_mid = view_2d();
283  p_mid = view_2d();
284  p_int = view_2d();
285  dens = view_2d();
286  omega = view_2d();
287  qv = view_2d();
288  qc = view_2d();
289  qc_copy = view_2d();
290  z_mid = view_2d();
291  z_int = view_2d();
292  shoc_s = view_2d();
293  tke = view_2d();
294  tke_copy = view_2d();
295  rrho = view_2d();
296  rrho_i = view_2d();
297  thv = view_2d();
298  dz = view_2d();
299  dse = view_2d();
300  zt_grid = view_2d();
301  zi_grid = view_2d();
302  wm_zt = view_2d();
303  inv_exner = view_2d();
304  thlm = view_2d();
305  qw = view_2d();
306  cloud_frac = view_2d();
307  cldfrac_liq = view_2d();
309 
311 
312  // SHOCPostprocess data structures
313  //=======================================================
314  shoc_ql2 = view_2d();
316 
317  // SHOC InputOutput data structures
318  //=======================================================
320  tk = view_2d();
321 
322  horiz_wind = view_3d();
323 
324  // SHOC Output data structures
325  //=======================================================
326  pblh = view_1d();
327  ustar = view_1d();
328  obklen = view_1d();
329 
330  tkh = view_2d();
331 
332  // SHOC HistoryOutput data structures
333  //=======================================================
334  shoc_mix = view_2d();
335  w_sec = view_2d();
336  thl_sec = view_2d();
337  qw_sec = view_2d();
338  qwthl_sec = view_2d();
339  wthl_sec = view_2d();
340  wqw_sec = view_2d();
341  wtke_sec = view_2d();
342  uw_sec = view_2d();
343  vw_sec = view_2d();
344  w3 = view_2d();
345  wqls_sec = view_2d();
346  brunt = view_2d();
347  isotropy = view_2d();
348  shoc_cond = view_2d();
349  shoc_evap = view_2d();
350 
351  // SHOC Miscellaneous data structures
352  //=======================================================
353  // Small kernels
354  se_b = view_1d();
355  ke_b = view_1d();
356  wv_b = view_1d();
357  wl_b = view_1d();
358  se_a = view_1d();
359  ke_a = view_1d();
360  wv_a = view_1d();
361  wl_a = view_1d();
362  kbfs = view_1d();
363  ustar2 = view_1d();
364  wstar = view_1d();
366 
367  rho_zt = view_2d();
368  shoc_qv = view_2d();
369  tabs = view_2d();
370  dz_zt = view_2d();
371  dz_zi = view_2d();
372 }

◆ finalize_impl()

void SHOCInterface::finalize_impl ( )
723 {
724  // Do nothing (per SHOCMacrophysics::finalize_impl())
725 
726  // Fill the AMReX MFs from Kokkos Views
728 
729  // Deallocate the buffer arrays
730  dealloc_buffers();
731 }
void kokkos_buffers_to_mf()
Definition: ERF_ShocInterface.cpp:533
void dealloc_buffers()
Definition: ERF_ShocInterface.cpp:268

◆ init_buffers()

void SHOCInterface::init_buffers ( )
protected

◆ initialize_impl()

void SHOCInterface::initialize_impl ( )
539 {
540  // For now, set z_int(i,nlevs) = z_surf = 0
541  const Real z_surf = 0.0;
542 
543  // Set preprocess variables
549 
550  // Input Variables:
551  input.zt_grid = zt_grid;
552  input.zi_grid = zi_grid;
553  input.pres = p_mid;
554  input.presi = p_int;
555  input.pdel = dens;
556  input.thv = thv;
557  input.w_field = wm_zt;
558  input.wthl_sfc = wpthlp_sfc;
559  input.wqw_sfc = wprtp_sfc;
560  input.uw_sfc = upwp_sfc;
561  input.vw_sfc = vpwp_sfc;
562  input.wtracer_sfc = wtracer_sfc;
563  input.inv_exner = inv_exner;
564  input.phis = phis;
565 
566  // Input/Output Variables
567  input_output.host_dse = shoc_s;
568  input_output.tke = tke_copy;
569  input_output.thetal = thlm;
570  input_output.qw = qw;
571  input_output.horiz_wind = horiz_wind;
572  input_output.wthv_sec = sgs_buoy_flux;
573  input_output.qtracers = qtracers;
574  input_output.tk = tk;
575  input_output.shoc_cldfrac = cldfrac_liq;
576  input_output.shoc_ql = qc_copy;
577 
578  // Output Variables
579  output.pblh = pblh;
580  output.shoc_ql2 = shoc_ql2;
581  output.tkh = tkh;
582  output.ustar = ustar;
583  output.obklen = obklen;
584 
585  // Output (diagnostic)
586  history_output.shoc_mix = shoc_mix;
587  history_output.isotropy = isotropy;
588  history_output.shoc_cond = shoc_cond;
589  history_output.shoc_evap = shoc_evap;
590  history_output.w_sec = w_sec;
591  history_output.thl_sec = thl_sec;
592  history_output.qw_sec = qw_sec;
593  history_output.qwthl_sec = qwthl_sec;
594  history_output.wthl_sec = wthl_sec;
595  history_output.wqw_sec = wqw_sec;
596  history_output.wtke_sec = wtke_sec;
597  history_output.uw_sec = uw_sec;
598  history_output.vw_sec = vw_sec;
599  history_output.w3 = w3;
600  history_output.wqls_sec = wqls_sec;
601  history_output.brunt = brunt;
602 
603  // Shoc small kernels
604  temporaries.se_b = se_b;
605  temporaries.ke_b = ke_b;
606  temporaries.wv_b = wv_b;
607  temporaries.wl_b = wl_b;
608  temporaries.se_a = se_a;
609  temporaries.ke_a = ke_a;
610  temporaries.wv_a = wv_a;
611  temporaries.wl_a = wl_a;
612  temporaries.kbfs = kbfs;
613  temporaries.ustar2 = ustar2;
614  temporaries.wstar = wstar;
615  temporaries.rho_zt = rho_zt;
616  temporaries.shoc_qv = shoc_qv;
617  temporaries.tabs = tabs;
618  temporaries.dz_zt = dz_zt;
619  temporaries.dz_zi = dz_zi;
620 
621  // Set postprocess variables
623  rrho, qv, qw, qc, qc_copy, tke,
626  T_mid, dse, z_mid, phis);
627 
628  // Maximum number of levels in pbl from surface
629  const int ntop_shoc = 0;
630  const int nbot_shoc = m_num_layers;
631  view_1d pref_mid("pref_mid", m_num_layers);
632  Spack* s_mem = reinterpret_cast<Spack*>(pref_mid.data());
633  uview_1d<Spack> pref_mid_um(s_mem, m_num_layers);
634  Kokkos::parallel_for(Kokkos::RangePolicy<>(0, m_num_layers), KOKKOS_LAMBDA (const int ilev)
635  {
636  pref_mid_um(ilev) = p_mid(0,ilev);
637  });
638  Kokkos::fence();
639  m_npbl = SHF::shoc_init(nbot_shoc,ntop_shoc,pref_mid_um);
640 
641  // Cell length for input dx and dy
642  view_1d cell_length_x("cell_length_x", m_num_cols);
643  view_1d cell_length_y("cell_length_y", m_num_cols);
644  Kokkos::deep_copy(cell_length_x, m_geom.CellSize(0));
645  Kokkos::deep_copy(cell_length_y, m_geom.CellSize(1));
646  input.dx = cell_length_x;
647  input.dy = cell_length_y;
648 }
SHOCPreprocess shoc_preprocess
Definition: ERF_ShocInterface.H:579
amrex::Geometry m_geom
Definition: ERF_ShocInterface.H:551
SHF::SHOCOutput output
Definition: ERF_ShocInterface.H:571
Int m_npbl
Definition: ERF_ShocInterface.H:538
SHF::SHOCInput input
Definition: ERF_ShocInterface.H:569
SHOCPostprocess shoc_postprocess
Definition: ERF_ShocInterface.H:580
typename SHF::Spack Spack
Definition: ERF_ShocInterface.H:43
SHF::SHOCInputOutput input_output
Definition: ERF_ShocInterface.H:570
SHF::SHOCHistoryOutput history_output
Definition: ERF_ShocInterface.H:572
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:387
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:252

◆ kokkos_buffers_to_mf()

void SHOCInterface::kokkos_buffers_to_mf ( )
534 {}

◆ mf_to_kokkos_buffers()

void SHOCInterface::mf_to_kokkos_buffers ( )
377 {
378  // Expose for device (shoc preprocess)
379  //=======================================================
380  auto T_mid_d = T_mid;
381  auto p_mid_d = p_mid;
382  auto p_int_d = p_int;
383  auto dens_d = dens;
384  auto omega_d = omega;
385  auto phis_d = phis;
386  auto surf_sens_flux_d = surf_sens_flux;
387  auto surf_evap_d = surf_evap;
388  auto surf_mom_flux_d = surf_mom_flux;
389  auto qtracers_d = qtracers;
390  auto qv_d = qv;
391  auto qc_d = qc;
392  auto qc_copy_d = qc_copy;
393  auto z_mid_d = z_mid;
394  auto z_int_d = z_int;
395  auto shoc_s_d = shoc_s;
396  auto tke_d = tke;
397  auto tke_copy_d = tke_copy;
398  auto rrho_d = rrho;
399  auto rrho_i_d = rrho_i;
400  auto thv_d = thv;
401  auto dz_d = dz;
402  auto dse_d = dse;
403  auto zt_grid_d = zt_grid;
404  auto zi_grid_d = zi_grid;
405  auto wpthlp_sfc_d = wpthlp_sfc;
406  auto wprtp_sfc_d = wprtp_sfc;
407  auto upwp_sfc_d = upwp_sfc;
408  auto vpwp_sfc_d = vpwp_sfc;
409  auto wtracer_sfc_d = wtracer_sfc;
410  auto wm_zt_d = wm_zt;
411  auto inv_exner_d = inv_exner;
412  auto thlm_d = thlm;
413  auto qw_d = qw;
414  auto cloud_frac_d = cloud_frac;
415  auto cldfrac_liq_d = cldfrac_liq;
416  auto cldfrac_liq_prev_d = cldfrac_liq_prev;
417 
418  int ncol = m_num_cols;
419  int nlay = m_num_layers;
420  Real dz = m_geom.CellSize(2);
421  bool moist = (m_cons->nComp() > RhoQ1_comp);
422  bool ice = (m_cons->nComp() > RhoQ3_comp);
423  auto ProbLoArr = m_geom.ProbLoArray();
424  for (MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
425  const auto& vbx = mfi.validbox();
426  const int nx = vbx.length(0);
427  const int imin = vbx.smallEnd(0);
428  const int jmin = vbx.smallEnd(1);
429  const int offset = m_col_offsets[mfi.index()];
430  const Array4<const Real>& cons_arr = m_cons->const_array(mfi);
431  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
432  Array4<const Real>{};
433  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
434  {
435  // map [i,j,k] 0-based to [icol, ilay] 0-based
436  const int icol = (j-jmin)*nx + (i-imin) + offset;
437  const int ilay = k;
438 
439  // EOS input (at CC)
440  Real r = cons_arr(i,j,k,Rho_comp);
441  Real rt = cons_arr(i,j,k,RhoTheta_comp);
442  Real qv = (moist) ? cons_arr(i,j,k,RhoQ1_comp)/r : 0.0;
443  Real qc = (moist) ? cons_arr(i,j,k,RhoQ2_comp)/r : 0.0;
444  Real qi = (ice) ? cons_arr(i,j,k,RhoQ3_comp)/r : 0.0;
445 
446  // EOS avg to z-face
447  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
448  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
449  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
450  Real r_avg = 0.5 * (r + r_lo);
451  Real rt_avg = 0.5 * (rt + rt_lo);
452  Real qv_avg = 0.5 * (qv + qv_lo);
453 
454  // Fill view1d
455  phis_d(icol) = 0.0;
456  surf_sens_flux_d(icol) = 0.0;
457  surf_evap_d(icol) = 0.0;
458  surf_mom_flux_d(icol,0) = 0.0;
459  surf_mom_flux_d(icol,1) = 0.0;
460 
461  // Fill view2d at CC
462  dens_d(icol,ilay) = r;
463  T_mid_d(icol,ilay) = getTgivenRandRTh(r, rt, qv);
464  p_mid_d(icol,ilay) = getPgivenRTh(rt, qv);
465  qv_d(icol,ilay) = qv;
466  qc_d(icol,ilay) = qc;
467  qc_copy_d(icol,ilay) = qc;
468  tke_d(icol,ilay) = std::max(cons_arr(i,j,k,RhoKE_comp)/r, 0.0);
469  tke_copy_d(icol,ilay) = std::max(cons_arr(i,j,k,RhoKE_comp)/r, 0.0);
470 
471  dz_d(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
472  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
473  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
474  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
475  z_mid_d(icol,ilay) = (z_arr) ? 0.125 * ( z_arr(i ,j ,k+1) + z_arr(i ,j ,k)
476  + z_arr(i+1,j ,k+1) + z_arr(i+1,j ,k)
477  + z_arr(i ,j+1,k+1) + z_arr(i ,j+1,k)
478  + z_arr(i+1,j+1,k+1) + z_arr(i+1,j+1,k) ) :
479  ProbLoArr[2] + (k + 0.5) * dz;
480 
481  // Fill view2d at w-face
482  p_int_d(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
483  z_int_d(icol,ilay) = (z_arr) ? 0.25 * ( z_arr(i ,j ,k) + z_arr(i+1,j ,k)
484  + z_arr(i ,j+1,k) + z_arr(i+1,j+1,k) ) :
485  ProbLoArr[2] + (k) * dz;
486  if (ilay==(nlay-1)) {
487  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
488  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
489  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,0.0) : 0.0;
490  r_avg = 0.5 * (r + r_hi);
491  rt_avg = 0.5 * (rt + rt_hi);
492  qv_avg = 0.5 * (qv + qv_hi);
493  p_int_d(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
494  z_int_d(icol,ilay+1) = (z_arr) ? 0.25 * ( z_arr(i ,j ,k+1) + z_arr(i+1,j ,k+1)
495  + z_arr(i ,j+1,k+1) + z_arr(i+1,j+1,k+1) ) :
496  ProbLoArr[2] + (k+1) * dz;
497  }
498 
499  // Questionable guesses
500  zt_grid_d(icol,ilay) = z_mid_d(icol,ilay); // Our thermo grid is the height grid?
501  zi_grid_d(icol,ilay) = z_int_d(icol,ilay); // Heights of interface grid?
502  cloud_frac_d(icol,ilay) = ((qc+qi)>0.0) ? 1. : 0.;
503  cldfrac_liq_d(icol,ilay) = (qc>0.0) ? 1. : 0.;
504  cldfrac_liq_prev_d(icol,ilay) = (qc>0.0) ? 1. : 0.;
505 
506  // Don't know
507  wpthlp_sfc_d(icol) = 0.0;
508  wprtp_sfc_d(icol) = 0.0;
509  upwp_sfc_d(icol) = 0.0;
510  vpwp_sfc_d(icol) = 0.0;
511  wtracer_sfc_d(icol,0) = 0.0;
512 
513  omega_d(icol,ilay) = 0.0;
514  shoc_s_d(icol,ilay) = 0.0;
515  rrho_d(icol,ilay) = 0.0;
516  rrho_i_d(icol,ilay) = 0.0;
517  thv_d(icol,ilay) = 0.0;
518  wm_zt_d(icol,ilay) = 0.0;
519  inv_exner_d(icol,ilay) = 0.0;
520  thlm_d(icol,ilay) = 0.0;
521  qw_d(icol,ilay) = 0.0;
522  dse(icol,ilay) = 0.0;
523 
524 
525  // Fill view3d
526  qtracers_d(icol,ilay,0) = 0.0;
527  });
528  }
529 }
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 getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ3_comp
Definition: ERF_IndexDefines.H:44
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::MultiFab * m_z_phys
Definition: ERF_ShocInterface.H:560
amrex::MultiFab * m_cons
Definition: ERF_ShocInterface.H:557
amrex::Vector< int > m_col_offsets
Definition: ERF_ShocInterface.H:532
Here is the call graph for this function:

◆ name()

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

◆ requested_buffer_size_in_bytes()

size_t SHOCInterface::requested_buffer_size_in_bytes ( ) const
protected

◆ run_impl()

void SHOCInterface::run_impl ( const Real  dt)
653 {
654  using TPF = ekat::TeamPolicyFactory<KT::ExeSpace>;
655 
656  EKAT_REQUIRE_MSG (dt<=300,
657  "Error! SHOC is intended to run with a timestep no longer than 5 minutes.\n"
658  " Please, reduce timestep (perhaps increasing subcycling iterations).\n");
659 
660  const auto nlev_packs = ekat::npack<Spack>(m_num_layers);
661  const auto scan_policy = TPF::get_thread_range_parallel_scan_team_policy(m_num_cols, nlev_packs);
662  const auto default_policy = TPF::get_default_team_policy(m_num_cols, nlev_packs);
663 
664  // Preprocessing of SHOC inputs. Kernel contains a parallel_scan,
665  // so a special TeamPolicy is required.
666  Kokkos::parallel_for("shoc_preprocess",
667  scan_policy,
669  Kokkos::fence();
670 
672  Kokkos::deep_copy(wtracer_sfc, 0);
673 
674  if (apply_tms) {
676  }
677 
678  if (check_flux_state) {
680  }
681 
682  // For now set the host timestep to the shoc timestep. This forces
683  // number of SHOC timesteps (nadv) to be 1.
684  // TODO: input parameter?
685  hdtime = dt;
686  m_nadv = std::max(static_cast<int>(round(hdtime/dt)),1);
687 
688  // Reset internal WSM variables.
689  workspace_mgr.reset_internals();
690 
691  // Run shoc main
694 #ifdef SCREAM_SHOC_SMALL_KERNELS
695  , temporaries
696 #endif
697  );
698 
699  // Postprocessing of SHOC outputs
700  Kokkos::parallel_for("shoc_postprocess",
701  default_policy,
703  Kokkos::fence();
704 
705  // Extra SHOC output diagnostics
706  if (runtime_options.extra_diags) {
707  Kokkos::deep_copy(shoc_mix,history_output.shoc_mix);
708  Kokkos::deep_copy(brunt,history_output.brunt);
709  Kokkos::deep_copy(w3,history_output.w3);
710  Kokkos::deep_copy(isotropy,history_output.isotropy);
711  Kokkos::deep_copy(wthl_sec,history_output.wthl_sec);
712  Kokkos::deep_copy(wqw_sec,history_output.wqw_sec);
713  Kokkos::deep_copy(uw_sec,history_output.uw_sec);
714  Kokkos::deep_copy(vw_sec,history_output.vw_sec);
715  Kokkos::deep_copy(qw_sec,history_output.qw_sec);
716  Kokkos::deep_copy(thl_sec,history_output.thl_sec);
717  }
718 }
ekat::WorkspaceManager< Spack, KT::Device > workspace_mgr
Definition: ERF_ShocInterface.H:583
Int m_nadv
Definition: ERF_ShocInterface.H:539
void check_flux_state_consistency(const double dt)
Definition: ERF_ShocInterface.cpp:751
Int hdtime
Definition: ERF_ShocInterface.H:542
void apply_turbulent_mountain_stress()
Definition: ERF_ShocInterface.cpp:735
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:242

◆ set_computed_group_impl()

void SHOCInterface::set_computed_group_impl ( )
protected

◆ set_grids()

void SHOCInterface::set_grids ( int &  level,
const amrex::BoxArray &  ba,
amrex::Geometry &  geom,
amrex::MultiFab *  cons,
amrex::MultiFab *  z_phys 
)
117 {
118  // Set data members that may change
119  m_lev = level;
120  m_geom = geom;
121  m_cons = cons;
122  m_z_phys = z_phys;
123 
124  // Ensure the boxes span klo -> khi
125  int klo = geom.Domain().smallEnd(2);
126  int khi = geom.Domain().bigEnd(2);
127 
128  // Reset vector of offsets for columnar data
129  m_num_layers = geom.Domain().length(2);
130 
131  m_num_cols = 0;
132  m_col_offsets.clear();
133  m_col_offsets.resize(int(ba.size()));
134  for (amrex::MFIter mfi(*m_cons); mfi.isValid(); ++mfi) {
135  const auto& vbx = mfi.validbox();
136  AMREX_ALWAYS_ASSERT_WITH_MESSAGE((klo == vbx.smallEnd(2)) &&
137  (khi == vbx.bigEnd(2)),
138  "Vertical decomposition with shoc is not allowed.");
139  int nx = vbx.length(0);
140  int ny = vbx.length(1);
141  m_col_offsets[mfi.index()] = m_num_cols;
142  m_num_cols += nx * ny;
143  }
144 
145  // Allocate the buffer arrays
146  alloc_buffers();
147 
148  // Fill the KOKKOS Views from AMReX MFs
150 }
void alloc_buffers()
Definition: ERF_ShocInterface.cpp:154
int m_lev
Definition: ERF_ShocInterface.H:545
void mf_to_kokkos_buffers()
Definition: ERF_ShocInterface.cpp:376
@ cons
Definition: ERF_IndexDefines.H:140

Member Data Documentation

◆ apply_tms

bool SHOCInterface::apply_tms = false
protected

◆ brunt

view_2d SHOCInterface::brunt
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

◆ cloud_frac

view_2d SHOCInterface::cloud_frac
protected

◆ dens

view_2d SHOCInterface::dens
protected

◆ dse

view_2d SHOCInterface::dse
protected

◆ dz

view_2d SHOCInterface::dz
protected

◆ dz_zi

view_2d SHOCInterface::dz_zi
protected

◆ dz_zt

view_2d SHOCInterface::dz_zt
protected

◆ hdtime

Int SHOCInterface::hdtime
protected

◆ history_output

SHF::SHOCHistoryOutput SHOCInterface::history_output
protected

◆ horiz_wind

view_3d SHOCInterface::horiz_wind
protected

◆ input

SHF::SHOCInput SHOCInterface::input
protected

◆ input_output

SHF::SHOCInputOutput SHOCInterface::input_output
protected

◆ inv_exner

view_2d SHOCInterface::inv_exner
protected

◆ inv_qc_relvar

view_2d SHOCInterface::inv_qc_relvar
protected

◆ isotropy

view_2d SHOCInterface::isotropy
protected

◆ kbfs

view_1d SHOCInterface::kbfs
protected

◆ ke_a

view_1d SHOCInterface::ke_a
protected

◆ ke_b

view_1d SHOCInterface::ke_b
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_lev

int SHOCInterface::m_lev
protected

◆ m_nadv

Int SHOCInterface::m_nadv
protected

◆ m_npbl

Int SHOCInterface::m_npbl
protected

◆ m_num_cols

Int SHOCInterface::m_num_cols
protected

◆ m_num_layers

Int SHOCInterface::m_num_layers
protected

◆ m_num_tracers

Int SHOCInterface::m_num_tracers = 1
protected

◆ m_num_vel_comp

Int SHOCInterface::m_num_vel_comp = 2
protected

◆ m_step

int SHOCInterface::m_step
protected

◆ m_z_phys

amrex::MultiFab* SHOCInterface::m_z_phys = 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

◆ qc

view_2d SHOCInterface::qc
protected

◆ qc_copy

view_2d SHOCInterface::qc_copy
protected

◆ qtracers

view_3d_strided SHOCInterface::qtracers
protected

◆ qv

view_2d SHOCInterface::qv
protected

◆ qw

view_2d SHOCInterface::qw
protected

◆ qw_sec

view_2d SHOCInterface::qw_sec
protected

◆ qwthl_sec

view_2d SHOCInterface::qwthl_sec
protected

◆ rho_zt

view_2d SHOCInterface::rho_zt
protected

◆ rrho

view_2d SHOCInterface::rrho
protected

◆ rrho_i

view_2d SHOCInterface::rrho_i
protected

◆ runtime_options

SHF::SHOCRuntime SHOCInterface::runtime_options
protected

◆ se_a

view_1d SHOCInterface::se_a
protected

◆ se_b

view_1d SHOCInterface::se_b
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

◆ shoc_ql2

view_2d SHOCInterface::shoc_ql2
protected

◆ shoc_qv

view_2d SHOCInterface::shoc_qv
protected

◆ shoc_s

view_2d SHOCInterface::shoc_s
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

◆ tabs

view_2d SHOCInterface::tabs
protected

◆ thl_sec

view_2d SHOCInterface::thl_sec
protected

◆ thlm

view_2d SHOCInterface::thlm
protected

◆ thv

view_2d SHOCInterface::thv
protected

◆ tk

view_2d SHOCInterface::tk
protected

◆ tke

view_2d SHOCInterface::tke
protected

◆ tke_copy

view_2d SHOCInterface::tke_copy
protected

◆ tkh

view_2d SHOCInterface::tkh
protected

◆ upwp_sfc

view_1d SHOCInterface::upwp_sfc
protected

◆ ustar

view_1d SHOCInterface::ustar
protected

◆ ustar2

view_1d SHOCInterface::ustar2
protected

◆ uw_sec

view_2d SHOCInterface::uw_sec
protected

◆ vpwp_sfc

view_1d SHOCInterface::vpwp_sfc
protected

◆ vw_sec

view_2d SHOCInterface::vw_sec
protected

◆ w3

view_2d SHOCInterface::w3
protected

◆ w_sec

view_2d SHOCInterface::w_sec
protected

◆ wl_a

view_1d SHOCInterface::wl_a
protected

◆ wl_b

view_1d SHOCInterface::wl_b
protected

◆ wm_zt

view_2d SHOCInterface::wm_zt
protected

◆ workspace_mgr

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

◆ wprtp_sfc

view_1d SHOCInterface::wprtp_sfc
protected

◆ wpthlp_sfc

view_1d SHOCInterface::wpthlp_sfc
protected

◆ wqls_sec

view_2d SHOCInterface::wqls_sec
protected

◆ wqw_sec

view_2d SHOCInterface::wqw_sec
protected

◆ wstar

view_1d SHOCInterface::wstar
protected

◆ wthl_sec

view_2d SHOCInterface::wthl_sec
protected

◆ wtke_sec

view_2d SHOCInterface::wtke_sec
protected

◆ wtracer_sfc

view_2d SHOCInterface::wtracer_sfc
protected

◆ wv_a

view_1d SHOCInterface::wv_a
protected

◆ wv_b

view_1d SHOCInterface::wv_b
protected

◆ z_int

view_2d SHOCInterface::z_int
protected

◆ z_mid

view_2d SHOCInterface::z_mid
protected

◆ zi_grid

view_2d SHOCInterface::zi_grid
protected

◆ zt_grid

view_2d SHOCInterface::zt_grid
protected

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