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_in, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
 
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_in = 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
 
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
 

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  // Construct parser object for following reads
10  ParmParse pp("erf.shoc");
11 
12  // Parse runtime inputs at start up
13  pp.get("lambda_low" , runtime_options.lambda_low );
14  pp.get("lambda_high" , runtime_options.lambda_high );
15  pp.get("lambda_slope" , runtime_options.lambda_slope );
16  pp.get("lambda_thresh" , runtime_options.lambda_thresh );
17  pp.get("thl2tune" , runtime_options.thl2tune );
18  pp.get("qw2tune" , runtime_options.qw2tune );
19  pp.get("qwthl2tune" , runtime_options.qwthl2tune );
20  pp.get("w2tune" , runtime_options.w2tune );
21  pp.get("length_fac" , runtime_options.length_fac );
22  pp.get("c_diag_3rd_mom" , runtime_options.c_diag_3rd_mom);
23  pp.get("coeff_kh" , runtime_options.Ckh );
24  pp.get("coeff_km" , runtime_options.Ckm );
25  pp.get("shoc_1p5tke" , runtime_options.shoc_1p5tke );
26  pp.get("extra_shoc_diags", runtime_options.extra_diags );
27 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:230
SHF::SHOCRuntime runtime_options
Definition: ERF_ShocInterface.H:572
Here is the call graph for this function:

Member Function Documentation

◆ alloc_buffers()

void SHOCInterface::alloc_buffers ( )
107 {
108  // SHOCPreprocess data structures
109  //=======================================================
110  phis = view_1d("Phis" , m_num_cols);
111  wpthlp_sfc = view_1d("Wpthlp_sfc" , m_num_cols);
112  wprtp_sfc = view_1d("Wprtp_sfc" , m_num_cols);
113  upwp_sfc = view_1d("Upwp_sfc" , m_num_cols);
114  vpwp_sfc = view_1d("Vpwp_sfc" , m_num_cols);
115  surf_sens_flux = view_1d("Sfc sens flux" , m_num_cols);
116  surf_evap = view_1d("Sfc latent flux", m_num_cols);
117  surf_mom_flux = sview_2d("Sfc mom flux" , m_num_cols, m_num_vel_comp);
118  wtracer_sfc = view_2d("Wtracer_sfc" , m_num_cols, m_num_tracers );
119 
120  T_mid = view_2d("T_mid" , m_num_cols, m_num_layers );
121  p_mid = view_2d("P_mid" , m_num_cols, m_num_layers );
122  p_int = view_2d("P_int" , m_num_cols, m_num_layers+1);
123  dens = view_2d("Density" , m_num_cols, m_num_layers );
124  omega = view_2d("Omega" , m_num_cols, m_num_layers );
125  qv = view_2d("Qv" , m_num_cols, m_num_layers );
126  qc = view_2d("Qc" , m_num_cols, m_num_layers );
127  qc_copy = view_2d("Qc_copy" , m_num_cols, m_num_layers );
128  z_mid = view_2d("Z_mid" , m_num_cols, m_num_layers );
129  z_int = view_2d("Z_int" , m_num_cols, m_num_layers+1);
130  shoc_s = view_2d("Shoc_s" , m_num_cols, m_num_layers );
131  tke = view_2d("Tke" , m_num_cols, m_num_layers );
132  tke_copy = view_2d("Tke_copy" , m_num_cols, m_num_layers );
133  rrho = view_2d("Rrho" , m_num_cols, m_num_layers );
134  rrho_i = view_2d("Rrho_i" , m_num_cols, m_num_layers );
135  thv = view_2d("Thv" , m_num_cols, m_num_layers );
136  dz = view_2d("dz" , m_num_cols, m_num_layers );
137  dse = view_2d("dse" , m_num_cols, m_num_layers );
138  zt_grid = view_2d("Zt_grid" , m_num_cols, m_num_layers );
139  zi_grid = view_2d("Zi_grid" , m_num_cols, m_num_layers );
140  wm_zt = view_2d("wm_zt" , m_num_cols, m_num_layers );
141  inv_exner = view_2d("Inv_exner" , m_num_cols, m_num_layers );
142  thlm = view_2d("Thlm" , m_num_cols, m_num_layers );
143  qw = view_2d("Qw" , m_num_cols, m_num_layers );
144  cloud_frac = view_2d("Cld_frac" , m_num_cols, m_num_layers );
145  cldfrac_liq = view_2d("Cld_frac_liq" , m_num_cols, m_num_layers );
146  cldfrac_liq_prev = view_2d("cld_frac_liq_prev" , m_num_cols, m_num_layers );
147 
148  // NOTE: Use layoutright format
149  Kokkos::LayoutStride layout(m_num_cols , m_num_cols*m_num_layers, // stride for dim0
150  m_num_layers , m_num_layers, // stride for dim1
151  m_num_tracers, 1 // stride for dim2
152  );
153  qtracers = view_3d_strided("Qtracers" , layout);
154 
155 
156  // SHOCPostprocess data structures
157  //=======================================================
159  inv_qc_relvar = view_2d("inv_qc_relvar", m_num_cols, m_num_layers);
160 
161  // SHOC InputOutput data structures
162  //=======================================================
163  sgs_buoy_flux = view_2d("sgs_buoy_flux", m_num_cols, m_num_layers);
164  tk = view_2d("eddy_diff_mom", m_num_cols, m_num_layers);
165 
167 
168  // SHOC Output data structures
169  //=======================================================
170  pblh = view_1d("pbl_height", m_num_cols);
171  ustar = view_1d("ustar" , m_num_cols);
172  obklen = view_1d("obklen" , m_num_cols);
173 
174  tkh = view_2d("eddy_diff_heat", m_num_cols, m_num_layers);
175 
176  // SHOC HistoryOutput data structures
177  //=======================================================
178  shoc_mix = view_2d("shoc_mix" , m_num_cols, m_num_layers);
179  w_sec = view_2d("w_sec" , m_num_cols, m_num_layers);
180  thl_sec = view_2d("thl_sec" , m_num_cols, m_num_layers);
181  qw_sec = view_2d("qw_sec" , m_num_cols, m_num_layers);
182  qwthl_sec = view_2d("qwthl_sec", m_num_cols, m_num_layers);
183  wthl_sec = view_2d("wthl_sec" , m_num_cols, m_num_layers);
184  wqw_sec = view_2d("wqw_sec" , m_num_cols, m_num_layers);
185  wtke_sec = view_2d("wtke_sec" , m_num_cols, m_num_layers);
186  uw_sec = view_2d("uw_sec" , m_num_cols, m_num_layers);
187  vw_sec = view_2d("vw_sec" , m_num_cols, m_num_layers);
188  w3 = view_2d("w3" , m_num_cols, m_num_layers);
189  wqls_sec = view_2d("wqls_sec" , m_num_cols, m_num_layers);
190  brunt = view_2d("brunt" , m_num_cols, m_num_layers);
191  isotropy = view_2d("isotropy" , m_num_cols, m_num_layers);
192  shoc_cond = view_2d("shoc_cond", m_num_cols, m_num_layers);
193  shoc_evap = view_2d("shoc_evap", m_num_cols, m_num_layers);
194 
195  // SHOC Miscellaneous data structures
196  //=======================================================
197 }
view_2d p_int
Definition: ERF_ShocInterface.H:597
view_2d wqw_sec
Definition: ERF_ShocInterface.H:652
view_2d inv_exner
Definition: ERF_ShocInterface.H:617
view_2d wtke_sec
Definition: ERF_ShocInterface.H:653
view_2d thl_sec
Definition: ERF_ShocInterface.H:648
view_2d tkh
Definition: ERF_ShocInterface.H:642
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:621
Int m_num_layers
Definition: ERF_ShocInterface.H:536
view_2d vw_sec
Definition: ERF_ShocInterface.H:655
view_2d thv
Definition: ERF_ShocInterface.H:610
view_2d shoc_cond
Definition: ERF_ShocInterface.H:660
typename ekat::KokkosTypes< KokkosDefaultDevice >::template view_2d< Real > sview_2d
Definition: ERF_ShocInterface.H:48
view_1d vpwp_sfc
Definition: ERF_ShocInterface.H:590
view_2d brunt
Definition: ERF_ShocInterface.H:658
Int m_num_vel_comp
Definition: ERF_ShocInterface.H:540
view_2d dens
Definition: ERF_ShocInterface.H:598
view_1d surf_sens_flux
Definition: ERF_ShocInterface.H:591
view_1d pblh
Definition: ERF_ShocInterface.H:639
Int m_num_cols
Definition: ERF_ShocInterface.H:535
view_2d sgs_buoy_flux
Definition: ERF_ShocInterface.H:633
view_2d p_mid
Definition: ERF_ShocInterface.H:596
view_2d uw_sec
Definition: ERF_ShocInterface.H:654
view_2d dse
Definition: ERF_ShocInterface.H:612
view_2d qc
Definition: ERF_ShocInterface.H:601
view_2d dz
Definition: ERF_ShocInterface.H:611
view_2d qwthl_sec
Definition: ERF_ShocInterface.H:650
view_2d shoc_mix
Definition: ERF_ShocInterface.H:646
view_2d qw_sec
Definition: ERF_ShocInterface.H:649
view_2d cldfrac_liq_prev
Definition: ERF_ShocInterface.H:622
view_2d qw
Definition: ERF_ShocInterface.H:619
view_2d cloud_frac
Definition: ERF_ShocInterface.H:620
view_2d omega
Definition: ERF_ShocInterface.H:599
typename SHF::view_3d_strided< Spack > view_3d_strided
Definition: ERF_ShocInterface.H:52
view_2d rrho_i
Definition: ERF_ShocInterface.H:609
view_1d obklen
Definition: ERF_ShocInterface.H:641
view_2d z_mid
Definition: ERF_ShocInterface.H:603
view_2d qc_copy
Definition: ERF_ShocInterface.H:602
view_2d rrho
Definition: ERF_ShocInterface.H:608
view_2d inv_qc_relvar
Definition: ERF_ShocInterface.H:629
view_1d surf_evap
Definition: ERF_ShocInterface.H:592
view_2d w_sec
Definition: ERF_ShocInterface.H:647
view_1d phis
Definition: ERF_ShocInterface.H:586
view_2d shoc_evap
Definition: ERF_ShocInterface.H:661
view_3d horiz_wind
Definition: ERF_ShocInterface.H:635
typename SHF::view_1d< Real > view_1d
Definition: ERF_ShocInterface.H:44
view_2d wqls_sec
Definition: ERF_ShocInterface.H:657
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:624
view_1d upwp_sfc
Definition: ERF_ShocInterface.H:589
view_1d ustar
Definition: ERF_ShocInterface.H:640
view_2d tke_copy
Definition: ERF_ShocInterface.H:607
view_2d wm_zt
Definition: ERF_ShocInterface.H:616
view_2d w3
Definition: ERF_ShocInterface.H:656
Int m_num_tracers
Definition: ERF_ShocInterface.H:539
sview_2d surf_mom_flux
Definition: ERF_ShocInterface.H:593
typename SHF::view_2d< SHF::Spack > view_2d
Definition: ERF_ShocInterface.H:46
view_2d zi_grid
Definition: ERF_ShocInterface.H:614
view_2d thlm
Definition: ERF_ShocInterface.H:618
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:615
view_2d shoc_ql2
Definition: ERF_ShocInterface.H:628
view_2d tk
Definition: ERF_ShocInterface.H:634
view_2d qv
Definition: ERF_ShocInterface.H:600
view_2d wthl_sec
Definition: ERF_ShocInterface.H:651
view_2d tke
Definition: ERF_ShocInterface.H:606
view_2d zt_grid
Definition: ERF_ShocInterface.H:613
view_2d T_mid
Definition: ERF_ShocInterface.H:595
view_2d shoc_s
Definition: ERF_ShocInterface.H:605
view_2d z_int
Definition: ERF_ShocInterface.H:604
typename SHF::view_3d< Spack > view_3d
Definition: ERF_ShocInterface.H:50
view_2d isotropy
Definition: ERF_ShocInterface.H:659
view_1d wpthlp_sfc
Definition: ERF_ShocInterface.H:587
view_1d wprtp_sfc
Definition: ERF_ShocInterface.H:588

◆ apply_turbulent_mountain_stress()

void SHOCInterface::apply_turbulent_mountain_stress ( )
protected

◆ check_flux_state_consistency()

void SHOCInterface::check_flux_state_consistency ( const double  dt)
protected

◆ dealloc_buffers()

void SHOCInterface::dealloc_buffers ( )
202 {
203  // SHOCPreprocess data structures
204  //=======================================================
205  phis = view_1d();
206  wpthlp_sfc = view_1d();
207  wprtp_sfc = view_1d();
208  upwp_sfc = view_1d();
209  vpwp_sfc = view_1d();
211  surf_evap = view_1d();
213  wtracer_sfc = view_2d();
214 
215  T_mid = view_2d();
216  p_mid = view_2d();
217  p_int = view_2d();
218  dens = view_2d();
219  omega = view_2d();
220  qv = view_2d();
221  qc = view_2d();
222  qc_copy = view_2d();
223  z_mid = view_2d();
224  z_int = view_2d();
225  shoc_s = view_2d();
226  tke = view_2d();
227  tke_copy = view_2d();
228  rrho = view_2d();
229  rrho_i = view_2d();
230  thv = view_2d();
231  dz = view_2d();
232  dse = view_2d();
233  zt_grid = view_2d();
234  zi_grid = view_2d();
235  wm_zt = view_2d();
236  inv_exner = view_2d();
237  thlm = view_2d();
238  qw = view_2d();
239  cloud_frac = view_2d();
240  cldfrac_liq = view_2d();
242 
244 
245  // SHOCPostprocess data structures
246  //=======================================================
247  shoc_ql2 = view_2d();
249 
250  // SHOC InputOutput data structures
251  //=======================================================
253  tk = view_2d();
254 
255  horiz_wind = view_3d();
256 
257  // SHOC Output data structures
258  //=======================================================
259  pblh = view_1d();
260  ustar = view_1d();
261  obklen = view_1d();
262 
263  tkh = view_2d();
264 
265  // SHOC HistoryOutput data structures
266  //=======================================================
267  shoc_mix = view_2d();
268  w_sec = view_2d();
269  thl_sec = view_2d();
270  qw_sec = view_2d();
271  qwthl_sec = view_2d();
272  wthl_sec = view_2d();
273  wqw_sec = view_2d();
274  wtke_sec = view_2d();
275  uw_sec = view_2d();
276  vw_sec = view_2d();
277  w3 = view_2d();
278  wqls_sec = view_2d();
279  brunt = view_2d();
280  isotropy = view_2d();
281  shoc_cond = view_2d();
282  shoc_evap = view_2d();
283 
284  // SHOC Miscellaneous data structures
285  //=======================================================
286 }

◆ finalize_impl()

void SHOCInterface::finalize_impl ( )
864 {
865  // Do nothing (per SHOCMacrophysics::finalize_impl())
866 
867  // Fill the AMReX MFs from Kokkos Views
869 
870  // Deallocate the buffer arrays
871  dealloc_buffers();
872 }
void kokkos_buffers_to_mf()
Definition: ERF_ShocInterface.cpp:447
void dealloc_buffers()
Definition: ERF_ShocInterface.cpp:201

◆ init_buffers()

void SHOCInterface::init_buffers ( )
protected

◆ initialize_impl()

void SHOCInterface::initialize_impl ( )
453 {
454  // For now, set z_int(i,nlevs) = z_surf = 0
455  const Real z_surf = 0.0;
456 
457  // Set preprocess variables
463 
464  // Input Variables:
465  input.zt_grid = zt_grid;
466  input.zi_grid = zi_grid;
467  input.pres = p_mid;
468  input.presi = p_int;
469  input.pdel = dens;
470  input.thv = thv;
471  input.w_field = wm_zt;
472  input.wthl_sfc = wpthlp_sfc;
473  input.wqw_sfc = wprtp_sfc;
474  input.uw_sfc = upwp_sfc;
475  input.vw_sfc = vpwp_sfc;
476  input.wtracer_sfc = wtracer_sfc;
477  input.inv_exner = inv_exner;
478  input.phis = phis;
479 
480  // Input/Output Variables
481  input_output.host_dse = shoc_s;
482  input_output.tke = tke_copy;
483  input_output.thetal = thlm;
484  input_output.qw = qw;
485  input_output.horiz_wind = horiz_wind;
486  input_output.wthv_sec = sgs_buoy_flux;
487  input_output.qtracers = qtracers;
488  input_output.tk = tk;
489  input_output.shoc_cldfrac = cldfrac_liq;
490  input_output.shoc_ql = qc_copy;
491 
492  // Output Variables
493  output.pblh = pblh;
494  output.shoc_ql2 = shoc_ql2;
495  output.tkh = tkh;
496  output.ustar = ustar;
497  output.obklen = obklen;
498 
499  // Output (diagnostic)
500  history_output.shoc_mix = shoc_mix;
501  history_output.isotropy = isotropy;
502  history_output.shoc_cond = shoc_cond;
503  history_output.shoc_evap = shoc_evap;
504  history_output.w_sec = w_sec;
505  history_output.thl_sec = thl_sec;
506  history_output.qw_sec = qw_sec;
507  history_output.qwthl_sec = qwthl_sec;
508  history_output.wthl_sec = wthl_sec;
509  history_output.wqw_sec = wqw_sec;
510  history_output.wtke_sec = wtke_sec;
511  history_output.uw_sec = uw_sec;
512  history_output.vw_sec = vw_sec;
513  history_output.w3 = w3;
514  history_output.wqls_sec = wqls_sec;
515  history_output.brunt = brunt;
516 
517  // Set postprocess variables
519  rrho, qv, qw, qc, qc_copy, tke,
522  T_mid, dse, z_mid, phis);
523 
524 #if 0
525  // Gather runtime options
526  runtime_options.lambda_low = m_params.get<double>("lambda_low");
527  runtime_options.lambda_high = m_params.get<double>("lambda_high");
528  runtime_options.lambda_slope = m_params.get<double>("lambda_slope");
529  runtime_options.lambda_thresh = m_params.get<double>("lambda_thresh");
530  runtime_options.thl2tune = m_params.get<double>("thl2tune");
531  runtime_options.qw2tune = m_params.get<double>("qw2tune");
532  runtime_options.qwthl2tune = m_params.get<double>("qwthl2tune");
533  runtime_options.w2tune = m_params.get<double>("w2tune");
534  runtime_options.length_fac = m_params.get<double>("length_fac");
535  runtime_options.c_diag_3rd_mom = m_params.get<double>("c_diag_3rd_mom");
536  runtime_options.Ckh = m_params.get<double>("coeff_kh");
537  runtime_options.Ckm = m_params.get<double>("coeff_km");
538  runtime_options.shoc_1p5tke = m_params.get<bool>("shoc_1p5tke");
539  runtime_options.extra_diags = m_params.get<bool>("extra_shoc_diags");
540 
541  // Initialize all of the structures that are passed to shoc_main in run_impl.
542  // Note: Some variables in the structures are not stored in the field manager. For these
543  // variables a local view is constructed.
544 
545  const auto& T_mid = get_field_out("T_mid").get_view<Spack**>();
546  const auto& p_mid = get_field_in("p_mid").get_view<const Spack**>();
547  const auto& p_int = get_field_in("p_int").get_view<const Spack**>();
548  const auto& pseudo_density = get_field_in("pseudo_density").get_view<const Spack**>();
549  const auto& omega = get_field_in("omega").get_view<const Spack**>();
550  const auto& surf_sens_flux = get_field_in("surf_sens_flux").get_view<const Real*>();
551  const auto& surf_evap = get_field_in("surf_evap").get_view<const Real*>();
552  const auto& surf_mom_flux = get_field_in("surf_mom_flux").get_view<const Real**>();
553  const auto& qtracers = get_group_out("turbulence_advected_tracers").m_monolithic_field->get_strided_view<Spack***>();
554  const auto& qc = get_field_out("qc").get_view<Spack**>();
555  const auto& qv = get_field_out("qv").get_view<Spack**>();
556  const auto& tke = get_field_out("tke").get_view<Spack**>();
557  const auto& cldfrac_liq = get_field_out("cldfrac_liq").get_view<Spack**>();
558  const auto& cldfrac_liq_prev = get_field_out("cldfrac_liq_prev").get_view<Spack**>();
559  const auto& sgs_buoy_flux = get_field_out("sgs_buoy_flux").get_view<Spack**>();
560  const auto& tk = get_field_out("eddy_diff_mom").get_view<Spack**>();
561  const auto& inv_qc_relvar = get_field_out("inv_qc_relvar").get_view<Spack**>();
562  const auto& phis = get_field_in("phis").get_view<const Real*>();
563 
564  // Alias local variables from temporary buffer
565  auto z_mid = m_buffer.z_mid;
566  auto z_int = m_buffer.z_int;
569  auto upwp_sfc = m_buffer.upwp_sfc;
570  auto vpwp_sfc = m_buffer.vpwp_sfc;
571  auto rrho = m_buffer.rrho;
572  auto rrho_i = m_buffer.rrho_i;
573  auto thv = m_buffer.thv;
574  auto dz = m_buffer.dz;
575  auto zt_grid = m_buffer.zt_grid;
576  auto zi_grid = m_buffer.zi_grid;
578  auto wm_zt = m_buffer.wm_zt;
580  auto thlm = m_buffer.thlm;
581  auto qw = m_buffer.qw;
582  auto dse = m_buffer.dse;
583  auto tke_copy = m_buffer.tke_copy;
584  auto qc_copy = m_buffer.qc_copy;
585  auto shoc_ql2 = m_buffer.shoc_ql2;
586 
587  // For now, set z_int(i,nlevs) = z_surf = 0
588  const Real z_surf = 0.0;
589 
590  // Some SHOC variables should be initialized uniformly if an Initial run
591  if (run_type==RunType::Initial){
592  Kokkos::deep_copy(sgs_buoy_flux,0.0);
593  Kokkos::deep_copy(tk,0.0);
594  Kokkos::deep_copy(tke,0.0004);
595  Kokkos::deep_copy(tke_copy,0.0004);
596  Kokkos::deep_copy(cldfrac_liq,0.0);
597  }
598 
604 
605  // Input Variables:
606  input.zt_grid = shoc_preprocess.zt_grid;
607  input.zi_grid = shoc_preprocess.zi_grid;
608  input.pres = p_mid;
609  input.presi = p_int;
610  input.pdel = pseudo_density;
611  input.thv = shoc_preprocess.thv;
612  input.w_field = shoc_preprocess.wm_zt;
613  input.wthl_sfc = shoc_preprocess.wpthlp_sfc;
614  input.wqw_sfc = shoc_preprocess.wprtp_sfc;
615  input.uw_sfc = shoc_preprocess.upwp_sfc;
616  input.vw_sfc = shoc_preprocess.vpwp_sfc;
617  input.wtracer_sfc = shoc_preprocess.wtracer_sfc;
618  input.inv_exner = shoc_preprocess.inv_exner;
619  input.phis = phis;
620 
621  // Input/Output Variables
626  input_output.horiz_wind = get_field_out("horiz_winds").get_view<Spack***>();
627  input_output.wthv_sec = sgs_buoy_flux;
629  input_output.tk = tk;
630  input_output.shoc_cldfrac = cldfrac_liq;
631  input_output.shoc_ql = qc_copy;
632 
633  // Output Variables
634  output.pblh = get_field_out("pbl_height").get_view<Real*>();
635  output.shoc_ql2 = shoc_ql2;
636  output.tkh = get_field_out("eddy_diff_heat").get_view<Spack**>();
637  output.ustar = get_field_out("ustar").get_view<Real*>();
638  output.obklen = get_field_out("obklen").get_view<Real*>();
639 
640  // Output (diagnostic)
641  history_output.shoc_mix = m_buffer.shoc_mix;
642  history_output.isotropy = m_buffer.isotropy;
643  if (m_params.get<bool>("extra_shoc_diags", false)) {
644  history_output.shoc_cond = get_field_out("shoc_cond").get_view<Spack**>();
645  history_output.shoc_evap = get_field_out("shoc_evap").get_view<Spack**>();
646  } else {
647  history_output.shoc_cond = m_buffer.unused;
648  history_output.shoc_evap = m_buffer.unused;
649  }
650  history_output.w_sec = get_field_out("w_variance").get_view<Spack**>();
651  history_output.thl_sec = m_buffer.thl_sec;
652  history_output.qw_sec = m_buffer.qw_sec;
653  history_output.qwthl_sec = m_buffer.qwthl_sec;
654  history_output.wthl_sec = m_buffer.wthl_sec;
655  history_output.wqw_sec = m_buffer.wqw_sec;
656  history_output.wtke_sec = m_buffer.wtke_sec;
657  history_output.uw_sec = m_buffer.uw_sec;
658  history_output.vw_sec = m_buffer.vw_sec;
660  history_output.wqls_sec = m_buffer.wqls_sec;
661  history_output.brunt = m_buffer.brunt;
662 
663 #ifdef SCREAM_SHOC_SMALL_KERNELS
664  temporaries.se_b = m_buffer.se_b;
665  temporaries.ke_b = m_buffer.ke_b;
666  temporaries.wv_b = m_buffer.wv_b;
667  temporaries.wl_b = m_buffer.wl_b;
668  temporaries.se_a = m_buffer.se_a;
669  temporaries.ke_a = m_buffer.ke_a;
670  temporaries.wv_a = m_buffer.wv_a;
671  temporaries.wl_a = m_buffer.wl_a;
672  temporaries.kbfs = m_buffer.kbfs;
673  temporaries.ustar2 = m_buffer.ustar2;
674  temporaries.wstar = m_buffer.wstar;
675 
676  temporaries.rho_zt = m_buffer.rho_zt;
677  temporaries.shoc_qv = m_buffer.shoc_qv;
678  temporaries.tabs = m_buffer.tabs;
679  temporaries.dz_zt = m_buffer.dz_zt;
680  temporaries.dz_zi = m_buffer.dz_zi;
681 #endif
682 
686  T_mid, dse, z_mid, phis);
687 
688  if (has_column_conservation_check()) {
689  const auto& vapor_flux = get_field_out("vapor_flux").get_view<Real*>();
690  const auto& water_flux = get_field_out("water_flux").get_view<Real*>();
691  const auto& ice_flux = get_field_out("ice_flux").get_view<Real*>();
692  const auto& heat_flux = get_field_out("heat_flux").get_view<Real*>();
694  vapor_flux, water_flux,
695  ice_flux, heat_flux);
696  }
697 
698  // Set field property checks for the fields in this process
699  using Interval = FieldWithinIntervalCheck;
700  using LowerBound = FieldLowerBoundCheck;
701  add_postcondition_check<Interval>(get_field_out("T_mid"),m_grid,100.0,500.0,false);
702  add_postcondition_check<Interval>(get_field_out("qc"),m_grid,0.0,0.1,false);
703  add_postcondition_check<Interval>(get_field_out("horiz_winds"),m_grid,-400.0,400.0,false);
704  add_postcondition_check<LowerBound>(get_field_out("pbl_height"),m_grid,0);
705  add_postcondition_check<Interval>(get_field_out("cldfrac_liq"),m_grid,0.0,1.0,false);
706  add_postcondition_check<LowerBound>(get_field_out("tke"),m_grid,0);
707  // For qv, ensure it doesn't get negative, by allowing repair of any neg value.
708  // TODO: use a repairable lb that clips only "small" negative values
709  add_postcondition_check<Interval>(get_field_out("qv"),m_grid,0,0.2,true);
710 
711  // Setup WSM for internal local variables
712  const auto nlev_packs = ekat::npack<Spack>(m_num_layers);
713  const auto nlevi_packs = ekat::npack<Spack>(m_num_layers+1);
714  const int n_wind_slots = ekat::npack<Spack>(2)*Spack::n;
715  const int n_trac_slots = ekat::npack<Spack>(m_num_tracers+3)*Spack::n;
716  const auto default_policy = ekat::ExeSpaceUtils<KT::ExeSpace>::get_default_team_policy(m_num_cols, nlev_packs);
717  workspace_mgr.setup(m_buffer.wsm_data, nlevi_packs, 14+(n_wind_slots+n_trac_slots), default_policy);
718 
719  // Calculate pref_mid, and use that to calculate
720  // maximum number of levels in pbl from surface
721  const auto pref_mid = m_buffer.pref_mid;
722  const auto s_pref_mid = ekat::scalarize(pref_mid);
723  const auto hyam = m_grid->get_geometry_data("hyam").get_view<const Real*>();
724  const auto hybm = m_grid->get_geometry_data("hybm").get_view<const Real*>();
725  const auto ps0 = C::P0;
726  const auto psref = ps0;
727  Kokkos::parallel_for(Kokkos::RangePolicy<>(0, m_num_layers), KOKKOS_LAMBDA (const int lev) {
728  s_pref_mid(lev) = ps0*hyam(lev) + psref*hybm(lev);
729  });
730  Kokkos::fence();
731 
732  const int ntop_shoc = 0;
733  const int nbot_shoc = m_num_layers;
734  m_npbl = SHF::shoc_init(nbot_shoc,ntop_shoc,pref_mid);
735 
736  // Compute cell length for input dx and dy.
737  const auto ncols = m_num_cols;
738  view_1d cell_length("cell_length", ncols);
739  if (m_grid->has_geometry_data("dx_short")) {
740  // In this case IOP is running with a planar geometry
741  auto dx = m_grid->get_geometry_data("dx_short").get_view<const Real,Host>()();
742  Kokkos::deep_copy(cell_length, dx*1000); // convert km -> m
743  } else {
744  const auto area = m_grid->get_geometry_data("area").get_view<const Real*>();
745  const auto lat = m_grid->get_geometry_data("lat").get_view<const Real*>();
746  Kokkos::parallel_for(ncols, KOKKOS_LAMBDA (const int icol) {
747  // For now, we are considering dy=dx. Here, we
748  // will need to compute dx/dy instead of cell_length
749  // if we have dy!=dx.
750  cell_length(icol) = PF::calculate_dx_from_area(area(icol),lat(icol));;
751  });
752  }
753  input.dx = cell_length;
754  input.dy = cell_length;
755 
756  // Initialize Kokkos SHOC pool allocator
757  // const size_t nvar = 300;
758  // const size_t nbnd = std::max(k_dist_sw_k->get_nband(),k_dist_sw_k->get_nband());
759  // const size_t ncol = gas_concs_k.ncol;
760  // const size_t nlay = gas_concs_k.nlay;
761  // auto my_size_ref = static_cast<unsigned long>(nvar * ncol * nlay * nbnd);
762 
763  // pool_t::init(my_size_ref);
764 
765  // We are now initialized!
766  initialized = true;
767 #endif
768 }
amrex::Real Real
Definition: ERF_ShocInterface.H:16
SHOCPreprocess shoc_preprocess
Definition: ERF_ShocInterface.H:578
SHF::SHOCOutput output
Definition: ERF_ShocInterface.H:570
ekat::WorkspaceManager< Spack, KT::Device > workspace_mgr
Definition: ERF_ShocInterface.H:582
Int m_npbl
Definition: ERF_ShocInterface.H:537
SHF::SHOCInput input
Definition: ERF_ShocInterface.H:568
SHOCPostprocess shoc_postprocess
Definition: ERF_ShocInterface.H:579
typename SHF::Spack Spack
Definition: ERF_ShocInterface.H:40
SHF::SHOCInputOutput input_output
Definition: ERF_ShocInterface.H:569
Buffer m_buffer
Definition: ERF_ShocInterface.H:565
SHF::SHOCHistoryOutput history_output
Definition: ERF_ShocInterface.H:571
bool initialized
Definition: ERF_RRTMGP_Interface.cpp:24
uview_2d< Spack > w3
Definition: ERF_ShocInterface.H:492
uview_2d< Spack > qc_copy
Definition: ERF_ShocInterface.H:479
uview_2d< Spack > dse
Definition: ERF_ShocInterface.H:477
uview_1d< Real > upwp_sfc
Definition: ERF_ShocInterface.H:446
uview_2d< Spack > vw_sec
Definition: ERF_ShocInterface.H:491
uview_1d< Real > wpthlp_sfc
Definition: ERF_ShocInterface.H:444
uview_2d< Spack > thv
Definition: ERF_ShocInterface.H:468
uview_2d< Spack > wqls_sec
Definition: ERF_ShocInterface.H:493
uview_2d< Spack > wthl_sec
Definition: ERF_ShocInterface.H:487
uview_2d< Spack > thlm
Definition: ERF_ShocInterface.H:475
uview_2d< Spack > uw_sec
Definition: ERF_ShocInterface.H:490
uview_2d< Spack > z_mid
Definition: ERF_ShocInterface.H:464
uview_1d< Real > vpwp_sfc
Definition: ERF_ShocInterface.H:447
uview_2d< Spack > rrho
Definition: ERF_ShocInterface.H:466
uview_2d< Spack > wtke_sec
Definition: ERF_ShocInterface.H:489
uview_2d< Spack > rrho_i
Definition: ERF_ShocInterface.H:467
uview_2d< Spack > unused
Definition: ERF_ShocInterface.H:463
Spack * wsm_data
Definition: ERF_ShocInterface.H:504
uview_2d< Spack > wtracer_sfc
Definition: ERF_ShocInterface.H:472
uview_2d< Spack > inv_exner
Definition: ERF_ShocInterface.H:474
uview_2d< Spack > qwthl_sec
Definition: ERF_ShocInterface.H:486
uview_2d< Spack > qw
Definition: ERF_ShocInterface.H:476
uview_2d< Spack > brunt
Definition: ERF_ShocInterface.H:494
uview_2d< Spack > shoc_mix
Definition: ERF_ShocInterface.H:481
uview_2d< Spack > tke_copy
Definition: ERF_ShocInterface.H:478
uview_2d< Spack > isotropy
Definition: ERF_ShocInterface.H:482
uview_2d< Spack > zt_grid
Definition: ERF_ShocInterface.H:470
uview_2d< Spack > thl_sec
Definition: ERF_ShocInterface.H:484
uview_2d< Spack > wm_zt
Definition: ERF_ShocInterface.H:473
uview_2d< Spack > qw_sec
Definition: ERF_ShocInterface.H:485
uview_2d< Spack > wqw_sec
Definition: ERF_ShocInterface.H:488
uview_2d< Spack > z_int
Definition: ERF_ShocInterface.H:465
uview_2d< Spack > shoc_ql2
Definition: ERF_ShocInterface.H:480
uview_1d< Real > wprtp_sfc
Definition: ERF_ShocInterface.H:445
uview_2d< Spack > zi_grid
Definition: ERF_ShocInterface.H:471
uview_1d< Spack > pref_mid
Definition: ERF_ShocInterface.H:461
uview_2d< Spack > dz
Definition: ERF_ShocInterface.H:469
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:386
void set_mass_and_energy_fluxes(const view_1d_const &surf_evap_, const view_1d_const &surf_sens_flux_, const view_1d &vapor_flux_, const view_1d &water_flux_, const view_1d &ice_flux_, const view_1d &heat_flux_)
Definition: ERF_ShocInterface.H:412
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:222
view_2d thlm
Definition: ERF_ShocInterface.H:244
view_2d zt_grid
Definition: ERF_ShocInterface.H:235
view_2d inv_exner
Definition: ERF_ShocInterface.H:243
view_2d wm_zt
Definition: ERF_ShocInterface.H:242
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:241
view_1d wprtp_sfc
Definition: ERF_ShocInterface.H:238
view_2d qw
Definition: ERF_ShocInterface.H:245
view_1d upwp_sfc
Definition: ERF_ShocInterface.H:239
view_1d wpthlp_sfc
Definition: ERF_ShocInterface.H:237
view_1d vpwp_sfc
Definition: ERF_ShocInterface.H:240
view_2d shoc_s
Definition: ERF_ShocInterface.H:228
view_2d tke_copy
Definition: ERF_ShocInterface.H:230
view_2d zi_grid
Definition: ERF_ShocInterface.H:236
view_2d thv
Definition: ERF_ShocInterface.H:233
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:251

◆ kokkos_buffers_to_mf()

void SHOCInterface::kokkos_buffers_to_mf ( )
448 {}

◆ mf_to_kokkos_buffers()

void SHOCInterface::mf_to_kokkos_buffers ( )
291 {
292  // Expose for device (shoc preprocess)
293  //=======================================================
294  auto T_mid_d = T_mid;
295  auto p_mid_d = p_mid;
296  auto p_int_d = p_int;
297  auto dens_d = dens;
298  auto omega_d = omega;
299  auto phis_d = phis;
300  auto surf_sens_flux_d = surf_sens_flux;
301  auto surf_evap_d = surf_evap;
302  auto surf_mom_flux_d = surf_mom_flux;
303  auto qtracers_d = qtracers;
304  auto qv_d = qv;
305  auto qc_d = qc;
306  auto qc_copy_d = qc_copy;
307  auto z_mid_d = z_mid;
308  auto z_int_d = z_int;
309  auto shoc_s_d = shoc_s;
310  auto tke_d = tke;
311  auto tke_copy_d = tke_copy;
312  auto rrho_d = rrho;
313  auto rrho_i_d = rrho_i;
314  auto thv_d = thv;
315  auto dz_d = dz;
316  auto dse_d = dse;
317  auto zt_grid_d = zt_grid;
318  auto zi_grid_d = zi_grid;
319  auto wpthlp_sfc_d = wpthlp_sfc;
320  auto wprtp_sfc_d = wprtp_sfc;
321  auto upwp_sfc_d = upwp_sfc;
322  auto vpwp_sfc_d = vpwp_sfc;
323  auto wtracer_sfc_d = wtracer_sfc;
324  auto wm_zt_d = wm_zt;
325  auto inv_exner_d = inv_exner;
326  auto thlm_d = thlm;
327  auto qw_d = qw;
328  auto cloud_frac_d = cloud_frac;
329  auto cldfrac_liq_d = cldfrac_liq;
330  auto cldfrac_liq_prev_d = cldfrac_liq_prev;
331 
332  int ncol = m_num_cols;
333  int nlay = m_num_layers;
334  Real dz = m_geom.CellSize(2);
335  bool moist = (m_cons_in->nComp() > RhoQ1_comp);
336  bool ice = (m_cons_in->nComp() > RhoQ3_comp);
337  auto ProbLoArr = m_geom.ProbLoArray();
338  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
339  const auto& vbx = mfi.validbox();
340  const int nx = vbx.length(0);
341  const int imin = vbx.smallEnd(0);
342  const int jmin = vbx.smallEnd(1);
343  const int offset = m_col_offsets[mfi.index()];
344  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
345  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
346  Array4<const Real>{};
347  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
348  {
349  // map [i,j,k] 0-based to [icol, ilay] 0-based
350  const int icol = (j-jmin)*nx + (i-imin) + offset;
351  const int ilay = k;
352 
353  // EOS input (at CC)
354  Real r = cons_arr(i,j,k,Rho_comp);
355  Real rt = cons_arr(i,j,k,RhoTheta_comp);
356  Real qv = (moist) ? cons_arr(i,j,k,RhoQ1_comp)/r : 0.0;
357  Real qc = (moist) ? cons_arr(i,j,k,RhoQ2_comp)/r : 0.0;
358  Real qi = (ice) ? cons_arr(i,j,k,RhoQ3_comp)/r : 0.0;
359 
360  // EOS avg to z-face
361  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
362  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
363  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
364  Real r_avg = 0.5 * (r + r_lo);
365  Real rt_avg = 0.5 * (rt + rt_lo);
366  Real qv_avg = 0.5 * (qv + qv_lo);
367 
368  // Fill view1d
369  phis_d(icol) = 0.0;
370  surf_sens_flux_d(icol) = 0.0;
371  surf_evap_d(icol) = 0.0;
372  surf_mom_flux_d(icol,0) = 0.0;
373  surf_mom_flux_d(icol,1) = 0.0;
374 
375  // Fill view2d at CC
376  dens_d(icol,ilay) = r;
377  T_mid_d(icol,ilay) = getTgivenRandRTh(r, rt, qv);
378  p_mid_d(icol,ilay) = getPgivenRTh(rt, qv);
379  qv_d(icol,ilay) = qv;
380  qc_d(icol,ilay) = qc;
381  qc_copy_d(icol,ilay) = qc;
382  tke_d(icol,ilay) = std::max(cons_arr(i,j,k,RhoKE_comp)/r, 0.0);
383  tke_copy_d(icol,ilay) = std::max(cons_arr(i,j,k,RhoKE_comp)/r, 0.0);
384 
385  dz_d(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
386  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
387  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
388  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
389  z_mid_d(icol,ilay) = (z_arr) ? 0.125 * ( z_arr(i ,j ,k+1) + z_arr(i ,j ,k)
390  + z_arr(i+1,j ,k+1) + z_arr(i+1,j ,k)
391  + z_arr(i ,j+1,k+1) + z_arr(i ,j+1,k)
392  + z_arr(i+1,j+1,k+1) + z_arr(i+1,j+1,k) ) :
393  ProbLoArr[2] + (k + 0.5) * dz;
394 
395  // Fill view2d at w-face
396  p_int_d(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
397  z_int_d(icol,ilay) = (z_arr) ? 0.25 * ( z_arr(i ,j ,k) + z_arr(i+1,j ,k)
398  + z_arr(i ,j+1,k) + z_arr(i+1,j+1,k) ) :
399  ProbLoArr[2] + (k) * dz;
400  if (ilay==(nlay-1)) {
401  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
402  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
403  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,0.0) : 0.0;
404  r_avg = 0.5 * (r + r_hi);
405  rt_avg = 0.5 * (rt + rt_hi);
406  qv_avg = 0.5 * (qv + qv_hi);
407  p_int_d(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
408  z_int_d(icol,ilay+1) = (z_arr) ? 0.25 * ( z_arr(i ,j ,k+1) + z_arr(i+1,j ,k+1)
409  + z_arr(i ,j+1,k+1) + z_arr(i+1,j+1,k+1) ) :
410  ProbLoArr[2] + (k+1) * dz;
411  }
412 
413  // Questionable guesses
414  zt_grid_d(icol,ilay) = z_mid_d(icol,ilay); // Our thermo grid is the height grid?
415  zi_grid_d(icol,ilay) = z_int_d(icol,ilay); // Heights of interface grid?
416  cloud_frac_d(icol,ilay) = ((qc+qi)>0.0) ? 1. : 0.;
417  cldfrac_liq_d(icol,ilay) = (qc>0.0) ? 1. : 0.;
418  cldfrac_liq_prev_d(icol,ilay) = (qc>0.0) ? 1. : 0.;
419 
420  // Don't know
421  wpthlp_sfc_d(icol) = 0.0;
422  wprtp_sfc_d(icol) = 0.0;
423  upwp_sfc_d(icol) = 0.0;
424  vpwp_sfc_d(icol) = 0.0;
425  wtracer_sfc_d(icol,0) = 0.0;
426 
427  omega_d(icol,ilay) = 0.0;
428  shoc_s_d(icol,ilay) = 0.0;
429  rrho_d(icol,ilay) = 0.0;
430  rrho_i_d(icol,ilay) = 0.0;
431  thv_d(icol,ilay) = 0.0;
432  wm_zt_d(icol,ilay) = 0.0;
433  inv_exner_d(icol,ilay) = 0.0;
434  thlm_d(icol,ilay) = 0.0;
435  qw_d(icol,ilay) = 0.0;
436  dse(icol,ilay) = 0.0;
437 
438 
439  // Fill view3d
440  qtracers_d(icol,ilay,0) = 0.0;
441  });
442  }
443 }
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::Geometry m_geom
Definition: ERF_ShocInterface.H:550
amrex::MultiFab * m_z_phys
Definition: ERF_ShocInterface.H:559
amrex::Vector< int > m_col_offsets
Definition: ERF_ShocInterface.H:531
amrex::MultiFab * m_cons_in
Definition: ERF_ShocInterface.H:556
Here is the call graph for this function:

◆ name()

std::string SHOCInterface::name ( ) const
inline
106 { 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)
773 {
774  EKAT_REQUIRE_MSG (dt<=300,
775  "Error! SHOC is intended to run with a timestep no longer than 5 minutes.\n"
776  " Please, reduce timestep (perhaps increasing subcycling iterations).\n");
777 
778 #if 0
779  const auto nlev_packs = ekat::npack<Spack>(m_num_layers);
780  const auto scan_policy = ekat::ExeSpaceUtils<KT::ExeSpace>::get_thread_range_parallel_scan_team_policy(m_num_cols, nlev_packs);
781  const auto default_policy = ekat::ExeSpaceUtils<KT::ExeSpace>::get_default_team_policy(m_num_cols, nlev_packs);
782 
783  // Preprocessing of SHOC inputs. Kernel contains a parallel_scan,
784  // so a special TeamPolicy is required.
785  Kokkos::parallel_for("shoc_preprocess",
786  scan_policy,
788  Kokkos::fence();
789 
791  Kokkos::deep_copy(wtracer_sfc, 0);
792 
793  if (m_params.get<bool>("apply_tms", false)) {
795  }
796 
797  if (m_params.get<bool>("check_flux_state_consistency", false)) {
799  }
800 
801  // For now set the host timestep to the shoc timestep. This forces
802  // number of SHOC timesteps (nadv) to be 1.
803  // TODO: input parameter?
804  hdtime = dt;
805  m_nadv = std::max(static_cast<int>(round(hdtime/dt)),1);
806 
807  // Reset internal WSM variables.
808  workspace_mgr.reset_internals();
809 
810  // Run shoc main
813 #ifdef SCREAM_SHOC_SMALL_KERNELS
814  , temporaries
815 #endif
816  );
817 
818  // Postprocessing of SHOC outputs
819  Kokkos::parallel_for("shoc_postprocess",
820  default_policy,
822  Kokkos::fence();
823 
824  // Extra SHOC output diagnostics
825  if (m_params.get<bool>("extra_shoc_diags", false)) {
826 
827  const auto& shoc_mix = get_field_out("shoc_mix").get_view<Spack**>();
828  Kokkos::deep_copy(shoc_mix,history_output.shoc_mix);
829 
830  const auto& brunt = get_field_out("brunt").get_view<Spack**>();
831  Kokkos::deep_copy(brunt,history_output.brunt);
832 
833  const auto& w3 = get_field_out("w3").get_view<Spack**>();
834  Kokkos::deep_copy(w3,history_output.w3);
835 
836  const auto& isotropy = get_field_out("isotropy").get_view<Spack**>();
837  Kokkos::deep_copy(isotropy,history_output.isotropy);
838 
839  const auto& wthl_sec = get_field_out("wthl_sec").get_view<Spack**>();
840  Kokkos::deep_copy(wthl_sec,history_output.wthl_sec);
841 
842  const auto& wqw_sec = get_field_out("wqw_sec").get_view<Spack**>();
843  Kokkos::deep_copy(wqw_sec,history_output.wqw_sec);
844 
845  const auto& uw_sec = get_field_out("uw_sec").get_view<Spack**>();
846  Kokkos::deep_copy(uw_sec,history_output.uw_sec);
847 
848  const auto& vw_sec = get_field_out("vw_sec").get_view<Spack**>();
849  Kokkos::deep_copy(vw_sec,history_output.vw_sec);
850 
851  const auto& qw_sec = get_field_out("qw_sec").get_view<Spack**>();
852  Kokkos::deep_copy(qw_sec,history_output.qw_sec);
853 
854  const auto& thl_sec = get_field_out("thl_sec").get_view<Spack**>();
855  Kokkos::deep_copy(thl_sec,history_output.thl_sec);
856 
857  } // Extra SHOC output diagnostics
858 #endif
859 }
Int m_nadv
Definition: ERF_ShocInterface.H:538
void check_flux_state_consistency(const double dt)
Int hdtime
Definition: ERF_ShocInterface.H:541
void apply_turbulent_mountain_stress()

◆ 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_in,
amrex::MultiFab *  z_phys,
amrex::MultiFab *  lat,
amrex::MultiFab *  lon 
)
69 {
70  // Set data members that may change
71  m_lev = level;
72  m_geom = geom;
73  m_cons_in = cons_in;
74  m_z_phys = z_phys;
75 
76  // Ensure the boxes span klo -> khi
77  int klo = geom.Domain().smallEnd(2);
78  int khi = geom.Domain().bigEnd(2);
79 
80  // Reset vector of offsets for columnar data
81  m_num_layers = geom.Domain().length(2);
82 
83  m_num_cols = 0;
84  m_col_offsets.clear();
85  m_col_offsets.resize(int(ba.size()));
86  for (amrex::MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
87  const auto& vbx = mfi.validbox();
88  AMREX_ALWAYS_ASSERT_WITH_MESSAGE((klo == vbx.smallEnd(2)) &&
89  (khi == vbx.bigEnd(2)),
90  "Vertical decomposition with shoc is not allowed.");
91  int nx = vbx.length(0);
92  int ny = vbx.length(1);
93  m_col_offsets[mfi.index()] = m_num_cols;
94  m_num_cols += nx * ny;
95  }
96 
97  // Allocate the buffer arrays
98  alloc_buffers();
99 
100  // Fill the KOKKOS Views from AMReX MFs
102 }
void alloc_buffers()
Definition: ERF_ShocInterface.cpp:106
int m_lev
Definition: ERF_ShocInterface.H:544
void mf_to_kokkos_buffers()
Definition: ERF_ShocInterface.cpp:290

Member Data Documentation

◆ brunt

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

◆ 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

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

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

◆ rrho

view_2d SHOCInterface::rrho
protected

◆ rrho_i

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

◆ shoc_ql2

view_2d SHOCInterface::shoc_ql2
protected

◆ shoc_s

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

◆ 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

◆ 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

◆ 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

◆ wthl_sec

view_2d SHOCInterface::wthl_sec
protected

◆ wtke_sec

view_2d SHOCInterface::wtke_sec
protected

◆ wtracer_sfc

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