ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ShocInterface.H
Go to the documentation of this file.
1 #ifndef ERF_SHOCINTERFACE_H
2 #define ERF_SHOCINTERFACE_H
3 
4 /*
5  * SHOC model interface to ERF -- we use this interface to link to the SHOC code
6  * as cloned from the E3SM git repository: https://github.com/E3SM-Project/E3SM/components/eamxx/src/physics/shoc
7  */
8 
9 #include "physics/shoc/shoc_constants.hpp"
10 #include "physics/shoc/shoc_functions.hpp"
11 #include "share/util/eamxx_common_physics_functions.hpp"
12 
13 #include "ekat_subview_utils.hpp"
14 #include "ekat_parameter_list.hpp"
15 #include "ekat_assert.hpp"
16 #include "ekat_team_policy_utils.hpp"
17 #include "ekat_reduction_utils.hpp"
18 
20 using Int=int;
21 
22 #include <string>
23 
24 #include <Kokkos_Core.hpp>
25 #include <AMReX_ParmParse.H>
26 #include <AMReX_MultiFabUtil.H>
27 #include <ERF_DataStruct.H>
28 #include <ERF_Kokkos.H>
29 #include <ERF_IndexDefines.H>
30 #include <ERF_EOS.H>
31 
32 //#include "share/atm_process/atmosphere_process.hpp"
33 //#include "share/atm_process/ATMBufferManager.hpp"
34 
36 {
37  using SHF = scream::shoc::Functions<Real, KokkosDefaultDevice>;
38  using PF = scream::PhysicsFunctions<KokkosDefaultDevice>;
39  using C = scream::physics::Constants<Real>;
40  using KT = ekat::KokkosTypes<KokkosDefaultDevice>;
41  using SC = scream::shoc::Constants<Real>;
42 
43  using Spack = typename SHF::Spack;
44  using IntSmallPack = typename SHF::IntSmallPack;
45  using Smask = typename SHF::Smask;
46  using view_1d_int = typename KT::template view_1d<Int>;
47  using view_1d = typename SHF::view_1d<Real>;
48  using view_1d_const = typename SHF::view_1d<const Real>;
49  using view_2d = typename SHF::view_2d<SHF::Spack>;
50  using view_2d_const = typename SHF::view_2d<const Spack>;
51  using sview_2d = typename KT::template view_2d<Real>;
52  using sview_2d_const = typename KT::template view_2d<const Real>;
53  using view_3d = typename SHF::view_3d<Spack>;
54  using view_3d_const = typename SHF::view_3d<const Spack>;
55  using view_3d_strided = typename SHF::view_3d_strided<Spack>;
56 
57  using WSM = ekat::WorkspaceManager<Spack, KT::Device>;
58 
59  template<typename ScalarT>
60  using uview_1d = ekat::Unmanaged<typename KT::template view_1d<ScalarT>>;
61  template<typename ScalarT>
62  using uview_2d = ekat::Unmanaged<typename KT::template view_2d<ScalarT>>;
63 
64 public:
65 
66  // Constructor
67  SHOCInterface (const int& lev,
68  SolverChoice& sc);
69 
70  // Set the grid
71  void
72  set_grids (int& level,
73  const amrex::BoxArray& ba,
74  amrex::Geometry& geom,
75  amrex::MultiFab* cons,
76  amrex::MultiFab* xvel,
77  amrex::MultiFab* yvel,
78  amrex::MultiFab* zvel,
79  amrex::Real* w_subsid,
80  amrex::MultiFab* tau13,
81  amrex::MultiFab* tau23,
82  amrex::MultiFab* hfx3,
83  amrex::MultiFab* qfx3,
84  amrex::MultiFab* eddyDiffs,
85  amrex::MultiFab* z_phys);
86 
87  // Initialize the implementation
88  void
89  initialize_impl ();
90 
91  // Run the implementation
92  void
93  run_impl (const Real dt);
94 
95  // Finalize the implementation
96  void
97  finalize_impl ();
98 
99  // Allocate buffer space
100  void
101  alloc_buffers ();
102 
103  // De-allocate buffer space
104  void
105  dealloc_buffers ();
106 
107  // Fill KOKKOS Views from AMReX MultiFabs
108  void
110 
111  // Fill AMReX MultiFabs from KOKKOS Views
112  void
114 
115  // The name of the subcomponent
116  std::string
117  name () const { return "shoc"; }
118 
119 #ifndef KOKKOS_ENABLE_CUDA
120  // Cuda requires methods enclosing __device__ lambda's to be public
121 protected:
122 #endif
123 
124  /*--------------------------------------------------------------------------------------------*/
125  // Most individual processes have a pre-processing step that constructs needed variables from
126  // the set of fields stored in the field manager. A structure like this defines those operations,
127  // which can then be called during run_impl in the main .cpp code.
128  // Structure to handle the local generation of data needed by shoc_main in run_impl
129  struct SHOCPreprocess {
130  SHOCPreprocess () = default;
131 
132  KOKKOS_INLINE_FUNCTION
133  void operator()(const Kokkos::TeamPolicy<KT::ExeSpace>::member_type& team) const
134  {
135  const int i = team.league_rank();
136 
137  const Real zvir = C::ZVIR;
138  const Real cpair = C::Cpair;
139  const Real ggr = C::gravit;
140  const Real inv_ggr = 1/ggr;
141  const Real mintke = SC::mintke;
142 
143  const int nlev_packs = ekat::npack<Spack>(nlev);
144 
145  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&] (const Int& k)
146  {
147  cldfrac_liq_prev(i,k)=cldfrac_liq(i,k);
148 
149  // Inverse of Exner. In non-rel builds, assert that exner != 0 when in range before computing.
150  const Spack exner = PF::exner_function(p_mid(i,k));
151  const Smask nonzero = (exner != 0);
152  EKAT_KERNEL_ASSERT((nonzero || !(ekat::range<IntSmallPack>(k*Spack::n) < nlev)).all());
153  inv_exner(i,k).set(nonzero, 1/exner);
154 
155  tke(i,k) = ekat::max(mintke, tke(i,k));
156 
157  // Tracers are updated as a group. The tracers tke and qc act as separate inputs to shoc_main()
158  // and are therefore updated differently to the tracers group's monolithic field. Here, we make
159  // a copy if each of these tracers and pass to shoc_main() so that changes to the tracer group
160  // does not alter tke or qc values. Then during post processing, we copy back correct values of
161  // tke and qc to tracer group in postprocessing.
162  // TODO: remove *_copy views once SHOC can request a subset of tracers.
163  tke_copy(i,k) = tke(i,k);
164  qc_copy(i,k) = qc(i,k);
165 
166  qw(i,k) = qv(i,k) + qc(i,k);
167 
168  // Temperature
169  // NOTE: theta_v (thv) is intentionally different from one in HOMME
170  const auto theta_zt = PF::calculate_theta_from_T(T_mid(i,k),p_mid(i,k));
171  thlm(i,k) = PF::calculate_thetal_from_theta(theta_zt,T_mid(i,k),qc(i,k));
172  thv(i,k) = theta_zt*(1 + zvir*qv(i,k) - qc(i,k));
173 
174  // Vertical layer thickness
175  //dz(i,k) = PF::calculate_dz(pseudo_density(i,k), p_mid(i,k), T_mid(i,k), qv(i,k));
176 
177  rrho(i,k) = inv_ggr*(pseudo_density(i,k)/dz(i,k));
178 
179  wm_zt(i,k) = -1*omega(i,k)/(rrho(i,k)*ggr);
180  });
181  team.team_barrier();
182 
183  // Compute vertical layer heights
184  const auto dz_s = ekat::subview(dz, i);
185  const auto z_int_s = ekat::subview(z_int, i);
186  const auto z_mid_s = ekat::subview(z_mid, i);
187  PF::calculate_z_int(team,nlev,dz_s,z_surf,z_int_s);
188  team.team_barrier();
189  PF::calculate_z_mid(team,nlev,z_int_s,z_mid_s);
190  team.team_barrier();
191 
192  const int nlevi_v = nlev/Spack::n;
193  const int nlevi_p = nlev%Spack::n;
194  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&] (const Int& k)
195  {
196  zt_grid(i,k) = z_mid(i,k) - z_int(i, nlevi_v)[nlevi_p];
197  zi_grid(i,k) = z_int(i,k) - z_int(i, nlevi_v)[nlevi_p];
198 
199  // Dry static energy
200  shoc_s(i,k) = PF::calculate_dse(T_mid(i,k),z_mid(i,k),phis(i));
201 
202  if (k+1 == nlev_packs) zi_grid(i,nlevi_v)[nlevi_p] = 0;
203  });
204  team.team_barrier();
205 
206  const auto zt_grid_s = ekat::subview(zt_grid, i);
207  const auto zi_grid_s = ekat::subview(zi_grid, i);
208  const auto rrho_s = ekat::subview(rrho, i);
209  const auto rrho_i_s = ekat::subview(rrho_i, i);
210  SHF::linear_interp(team,zt_grid_s,zi_grid_s,rrho_s,rrho_i_s,nlev,nlev+1,0);
211  team.team_barrier();
212 
213  const auto exner_int = PF::exner_function(p_int(i,nlevi_v)[nlevi_p]);
214  const auto inv_exner_int_surf = 1/exner_int;
215 
216  wpthlp_sfc(i) = (surf_sens_flux(i)/(cpair*rrho_i(i,nlevi_v)[nlevi_p]))*inv_exner_int_surf;
217  wprtp_sfc(i) = surf_evap(i)/rrho_i(i,nlevi_v)[nlevi_p];
218  upwp_sfc(i) = surf_mom_flux(i,0)/rrho_i(i,nlevi_v)[nlevi_p];
219  vpwp_sfc(i) = surf_mom_flux(i,1)/rrho_i(i,nlevi_v)[nlevi_p];
220  } // operator
221 
222  // Local variables
223  int ncol, nlev;
261 
262  // Assigning local variables
263  void set_variables(const int ncol_, const int nlev_,
264  const Real z_surf_,
265  const view_2d_const& T_mid_, const view_2d_const& p_mid_, const view_2d_const& p_int_, const view_2d_const& pseudo_density_,
266  const view_2d_const& omega_,
267  const view_1d_const& phis_, const view_1d_const& surf_sens_flux_, const view_1d_const& surf_evap_,
268  const sview_2d_const& surf_mom_flux_,
269  const view_3d_strided& qtracers_,
270  const view_2d& qv_, const view_2d_const& qc_, const view_2d& qc_copy_,
271  const view_2d& tke_, const view_2d& tke_copy_,
272  const view_2d& z_mid_, const view_2d& z_int_,
273  const view_2d& dse_, const view_2d& rrho_, const view_2d& rrho_i_,
274  const view_2d& thv_, const view_2d& dz_,const view_2d& zt_grid_,const view_2d& zi_grid_, const view_1d& wpthlp_sfc_,
275  const view_1d& wprtp_sfc_,const view_1d& upwp_sfc_,const view_1d& vpwp_sfc_, const view_2d& wtracer_sfc_,
276  const view_2d& wm_zt_,const view_2d& inv_exner_,const view_2d& thlm_,const view_2d& qw_,
277  const view_2d& cldfrac_liq_, const view_2d& cldfrac_liq_prev_)
278  {
279  ncol = ncol_;
280  nlev = nlev_;
281  z_surf = z_surf_;
282  // IN
283  T_mid = T_mid_;
284  p_mid = p_mid_;
285  p_int = p_int_;
286  pseudo_density = pseudo_density_;
287  omega = omega_;
288  phis = phis_;
289  surf_sens_flux = surf_sens_flux_;
290  surf_evap = surf_evap_;
291  surf_mom_flux = surf_mom_flux_;
292  qv = qv_;
293  // OUT
294  qtracers = qtracers_;
295  qc = qc_;
296  qc_copy = qc_copy_;
297  shoc_s = dse_;
298  tke = tke_;
299  tke_copy = tke_copy_;
300  z_mid = z_mid_;
301  z_int = z_int_;
302  rrho = rrho_;
303  rrho_i = rrho_i_;
304  thv = thv_;
305  dz = dz_;
306  zt_grid = zt_grid_;
307  zi_grid = zi_grid_;
308  wpthlp_sfc = wpthlp_sfc_;
309  wprtp_sfc = wprtp_sfc_;
310  upwp_sfc = upwp_sfc_;
311  vpwp_sfc = vpwp_sfc_;
312  wtracer_sfc = wtracer_sfc_;
313  wm_zt = wm_zt_;
314  inv_exner = inv_exner_;
315  thlm = thlm_;
316  qw = qw_;
317  cldfrac_liq=cldfrac_liq_;
318  cldfrac_liq_prev=cldfrac_liq_prev_;
319  } // set_variables
320  }; // SHOCPreprocess
321  /* --------------------------------------------------------------------------------------------*/
322 
323  /*--------------------------------------------------------------------------------------------*/
324  // Structure to handle the generation of data needed by the rest of the model based on output from
325  // shoc_main.
327  SHOCPostprocess () = default;
328 
329  KOKKOS_INLINE_FUNCTION
330  void operator()(const Kokkos::TeamPolicy<KT::ExeSpace>::member_type& team) const
331  {
332  const int i = team.league_rank();
333 
334  const Real inv_qc_relvar_max = 10;
335  const Real inv_qc_relvar_min = 0.001;
336 
337  const int nlev_packs = ekat::npack<Spack>(nlev);
338  Kokkos::parallel_for(Kokkos::TeamVectorRange(team, nlev_packs), [&] (const Int& k)
339  {
340  // See comment in SHOCPreprocess::operator() about the necessity of *_copy views
341  tke(i,k) = tke_copy(i,k);
342  qc(i,k) = qc_copy(i,k);
343 
344  qv(i,k) = qw(i,k) - qc(i,k);
345 
346  cldfrac_liq(i,k) = ekat::min(cldfrac_liq(i,k), 1);
347 
348  //P3 uses inv_qc_relvar, P3 is using dry mmrs, but
349  //wet<->dry conversion is a constant factor that cancels out in mean(qc)^2/mean(qc'*qc').
350  inv_qc_relvar(i,k) = 1;
351  const auto condition = (qc(i,k) != 0 && qc2(i,k) != 0);
352  if (condition.any()) {
353  inv_qc_relvar(i,k).set(condition,
354  ekat::min(inv_qc_relvar_max,
355  ekat::max(inv_qc_relvar_min,
356  ekat::square(qc(i,k))/qc2(i,k))));
357  }
358 
359  // Temperature
360  const Spack dse_ik(dse(i,k));
361  const Spack z_mid_ik(z_mid(i,k));
362  const Real phis_i(phis(i));
363  T_mid(i,k) = PF::calculate_temperature_from_dse(dse_ik,z_mid_ik,phis_i);
364 
365  });
366 
367  // If necessary, set appropriate boundary fluxes for energy and mass conservation checks.
368  // Any boundary fluxes not included in SHOC interface are set to 0.
370  vapor_flux(i) = surf_evap(i);
371  water_flux(i) = 0.0;
372  ice_flux(i) = 0.0;
373  heat_flux(i) = surf_sens_flux(i);
374  }
375  } // operator
376 
377  // Local variables
378  int ncol, nlev;
396 
397  // Assigning local variables
398  void set_variables (const int ncol_, const int nlev_,
399  const view_2d_const& rrho_,
400  const view_2d& qv_, const view_2d_const& qw_, const view_2d& qc_, const view_2d_const& qc_copy_,
401  const view_2d& tke_, const view_2d_const& tke_copy_, const view_3d_strided& qtracers_, const view_2d_const& qc2_,
402  const view_2d& cldfrac_liq_, const view_2d& inv_qc_relvar_,
403  const view_2d& T_mid_, const view_2d_const& dse_, const view_2d_const& z_mid_, const view_1d_const phis_)
404  {
405  ncol = ncol_;
406  nlev = nlev_;
407  rrho = rrho_;
408  qv = qv_;
409  qw = qw_;
410  qc = qc_;
411  qc_copy = qc_copy_;
412  tke = tke_;
413  tke_copy = tke_copy_;
414  qtracers = qtracers_;
415  qc2 = qc2_;
416  cldfrac_liq = cldfrac_liq_;
417  inv_qc_relvar = inv_qc_relvar_;
418  T_mid = T_mid_;
419  dse = dse_;
420  z_mid = z_mid_;
421  phis = phis_;
422  } // set_variables
423 
424  void set_mass_and_energy_fluxes (const view_1d_const& surf_evap_, const view_1d_const& surf_sens_flux_,
425  const view_1d& vapor_flux_, const view_1d& water_flux_,
426  const view_1d& ice_flux_, const view_1d& heat_flux_)
427  {
429  surf_evap = surf_evap_;
430  surf_sens_flux = surf_sens_flux_;
431  vapor_flux = vapor_flux_;
432  water_flux = water_flux_;
433  ice_flux = ice_flux_;
434  heat_flux = heat_flux_;
435  }
436  }; // SHOCPostprocess
437  /* --------------------------------------------------------------------------------------------*/
438 
439  // Structure for storing local variables initialized using the ATMBufferManager
440  struct Buffer {
441 #ifndef SCREAM_SHOC_SMALL_KERNELS
442  static constexpr int num_1d_scalar_ncol = 4;
443 #else
444  static constexpr int num_1d_scalar_ncol = 15;
445 #endif
446  static constexpr int num_1d_scalar_nlev = 1;
447 #ifndef SCREAM_SHOC_SMALL_KERNELS
448  static constexpr int num_2d_vector_mid = 19;
449  static constexpr int num_2d_vector_int = 12;
450 #else
451  static constexpr int num_2d_vector_mid = 23;
452  static constexpr int num_2d_vector_int = 13;
453 #endif
454  static constexpr int num_2d_vector_tr = 1;
455 
460 #ifdef SCREAM_SHOC_SMALL_KERNELS
461  uview_1d<Real> se_b;
462  uview_1d<Real> ke_b;
463  uview_1d<Real> wv_b;
464  uview_1d<Real> wl_b;
465  uview_1d<Real> se_a;
466  uview_1d<Real> ke_a;
467  uview_1d<Real> wv_a;
468  uview_1d<Real> wl_a;
469  uview_1d<Real> kbfs;
470  uview_1d<Real> ustar2;
471  uview_1d<Real> wstar;
472 #endif
474 
475  uview_2d<Spack> unused; // Placeholder for unused views
507 #ifdef SCREAM_SHOC_SMALL_KERNELS
508  uview_2d<Spack> rho_zt;
509  uview_2d<Spack> shoc_qv;
511  uview_2d<Spack> dz_zt;
512  uview_2d<Spack> dz_zi;
514 #endif
515 
517  }; // BUFFER
518  /* --------------------------------------------------------------------------------------------*/
519 
520 #ifndef KOKKOS_ENABLE_CUDA
521  // Cuda requires methods enclosing __device__ lambda's to be public
522 protected:
523 #endif
524 
525  // Update flux (if necessary)
526  void
527  check_flux_state_consistency (const double dt);
528 
529  // Apply TMS drag coeff to shoc_main inputs (if necessary)
530  void
532 
533 protected:
534 
535  // SHOC updates the 'tracers' group.
536  void
538 
539  // Computes total number of bytes needed for local variables
540  size_t
542 
543  // Set local variables
544  void
545  init_buffers ();
546 
547  // Offsets for MultiFab <-> KOKKOS transfer
548  amrex::Vector<int> m_col_offsets;
549 
550  // Keep track of field dimensions and other scalar values
551  // needed in shoc_main
559 
560  // Grid level
561  int m_lev;
562 
563  // Step number
564  int m_step;
565 
566  // Geometry at given level
567  amrex::Geometry m_geom;
568 
569  // Boxarray at given level
570  amrex::BoxArray m_ba;
571 
572  // Pointer to the CC conserved vars
573  amrex::MultiFab* m_cons = nullptr;
574 
575  // Pointer to the FC velocity vars
576  amrex::MultiFab* m_xvel = nullptr;
577  amrex::MultiFab* m_yvel = nullptr;
578  amrex::MultiFab* m_zvel = nullptr;
579 
580  // Pointer to the w-subsidence velocity
582 
583  // Pointer to the surface stresses
584  amrex::MultiFab* m_tau13 = nullptr;
585  amrex::MultiFab* m_tau23 = nullptr;
586  amrex::MultiFab* m_hfx3 = nullptr;
587  amrex::MultiFab* m_qfx3 = nullptr;
588 
589  // Pointer to the diffusion coeffs
590  amrex::MultiFab* m_mu = nullptr;
591 
592  // Pointer to the terrain heights
593  amrex::MultiFab* m_z_phys = nullptr;
594 
595  // Flag for first step
596  bool m_first_step = true;
597 
598  // Struct which contains local variables
600 
601  // Store the structures for each argument to shoc_main;
602  SHF::SHOCInput input;
603  SHF::SHOCInputOutput input_output;
604  SHF::SHOCOutput output;
605  SHF::SHOCHistoryOutput history_output;
606  SHF::SHOCRuntime runtime_options;
607 #ifdef SCREAM_SHOC_SMALL_KERNELS
608  SHF::SHOCTemporaries temporaries;
609 #endif
610 
611  // Structures which compute pre/post process
614 
615  // WSM for internal local variables
616  ekat::WorkspaceManager<Spack, KT::Device> workspace_mgr;
617 
618  // Contiguous memory buffer view
619  //=======================================================
621 
622  // Runtime options
623  //=======================================================
624  bool apply_tms = true;
625  bool check_flux_state = false;
626  bool extra_shoc_diags = false;
628 
629  // Interface data structures
630  //=======================================================
638 
639  // Input data structures
640  //=======================================================
645 
646  // Input/Output data structures
647  //=======================================================
654 
655  // Output data structures
656  //=======================================================
664 
665  // Extra diagnostic data structures
666  //=======================================================
672 
680 
681  // Tracer data structures
682  //=======================================================
684 
685  // Boundary flux data structures
686  //=======================================================
691 }; // SHOCInterface
692 #endif
@ tau23
Definition: ERF_DataStruct.H:30
@ tau13
Definition: ERF_DataStruct.H:30
int Int
Definition: ERF_ShocInterface.H:20
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_ShocInterface.H:36
view_2d p_int
Definition: ERF_ShocInterface.H:642
typename SHF::view_1d< const Real > view_1d_const
Definition: ERF_ShocInterface.H:48
bool check_flux_state
Definition: ERF_ShocInterface.H:625
SHOCPreprocess shoc_preprocess
Definition: ERF_ShocInterface.H:612
view_2d wqw_sec
Definition: ERF_ShocInterface.H:675
scream::shoc::Constants< Real > SC
Definition: ERF_ShocInterface.H:41
view_1d surf_drag_coeff_tms
Definition: ERF_ShocInterface.H:637
amrex::Geometry m_geom
Definition: ERF_ShocInterface.H:567
void set_computed_group_impl()
void alloc_buffers()
Definition: ERF_ShocInterface.cpp:186
view_2d thl_sec
Definition: ERF_ShocInterface.H:674
view_2d tkh
Definition: ERF_ShocInterface.H:659
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:651
typename KT::template view_2d< Real > sview_2d
Definition: ERF_ShocInterface.H:51
amrex::MultiFab * m_z_phys
Definition: ERF_ShocInterface.H:593
amrex::MultiFab * m_qfx3
Definition: ERF_ShocInterface.H:587
void init_buffers()
Definition: ERF_ShocInterface.cpp:617
Int m_num_layers
Definition: ERF_ShocInterface.H:553
view_2d vw_sec
Definition: ERF_ShocInterface.H:678
amrex::MultiFab * m_yvel
Definition: ERF_ShocInterface.H:577
view_1d ice_flux
Definition: ERF_ShocInterface.H:689
size_t requested_buffer_size_in_bytes() const
Definition: ERF_ShocInterface.cpp:592
view_2d shoc_cond
Definition: ERF_ShocInterface.H:670
view_2d brunt
Definition: ERF_ShocInterface.H:667
bool extra_shoc_diags
Definition: ERF_ShocInterface.H:626
view_1d heat_flux
Definition: ERF_ShocInterface.H:690
Int m_num_vel_comp
Definition: ERF_ShocInterface.H:557
bool apply_tms
Definition: ERF_ShocInterface.H:624
view_1d surf_sens_flux
Definition: ERF_ShocInterface.H:632
view_1d pblh
Definition: ERF_ShocInterface.H:657
amrex::MultiFab * m_zvel
Definition: ERF_ShocInterface.H:578
void set_grids(int &level, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons, amrex::MultiFab *xvel, amrex::MultiFab *yvel, amrex::MultiFab *zvel, amrex::Real *w_subsid, amrex::MultiFab *tau13, amrex::MultiFab *tau23, amrex::MultiFab *hfx3, amrex::MultiFab *qfx3, amrex::MultiFab *eddyDiffs, amrex::MultiFab *z_phys)
Definition: ERF_ShocInterface.cpp:114
Int m_num_cols
Definition: ERF_ShocInterface.H:552
view_2d sgs_buoy_flux
Definition: ERF_ShocInterface.H:649
SHF::SHOCOutput output
Definition: ERF_ShocInterface.H:604
typename KT::template view_1d< Int > view_1d_int
Definition: ERF_ShocInterface.H:46
ekat::WorkspaceManager< Spack, KT::Device > workspace_mgr
Definition: ERF_ShocInterface.H:616
view_2d p_mid
Definition: ERF_ShocInterface.H:641
view_2d uw_sec
Definition: ERF_ShocInterface.H:677
Int m_nadv
Definition: ERF_ShocInterface.H:555
view_2d qc
Definition: ERF_ShocInterface.H:653
void run_impl(const Real dt)
Definition: ERF_ShocInterface.cpp:851
amrex::MultiFab * m_mu
Definition: ERF_ShocInterface.H:590
Int m_npbl
Definition: ERF_ShocInterface.H:554
amrex::MultiFab * m_tau23
Definition: ERF_ShocInterface.H:585
view_2d shoc_mix
Definition: ERF_ShocInterface.H:668
view_2d qw_sec
Definition: ERF_ShocInterface.H:676
view_2d cldfrac_liq_prev
Definition: ERF_ShocInterface.H:661
scream::physics::Constants< Real > C
Definition: ERF_ShocInterface.H:39
SHF::SHOCInput input
Definition: ERF_ShocInterface.H:602
SHOCInterface(const int &lev, SolverChoice &sc)
Definition: ERF_ShocInterface.cpp:6
view_2d omega
Definition: ERF_ShocInterface.H:631
typename SHF::view_2d< const Spack > view_2d_const
Definition: ERF_ShocInterface.H:50
ekat::Unmanaged< typename KT::template view_1d< ScalarT > > uview_1d
Definition: ERF_ShocInterface.H:60
typename SHF::Smask Smask
Definition: ERF_ShocInterface.H:45
typename SHF::view_3d_strided< Spack > view_3d_strided
Definition: ERF_ShocInterface.H:55
view_1d vapor_flux
Definition: ERF_ShocInterface.H:687
view_1d obklen
Definition: ERF_ShocInterface.H:663
void kokkos_buffers_to_mf()
Definition: ERF_ShocInterface.cpp:485
scream::shoc::Functions< Real, KokkosDefaultDevice > SHF
Definition: ERF_ShocInterface.H:37
view_2d pseudo_dens
Definition: ERF_ShocInterface.H:643
SHF::SHOCRuntime runtime_options
Definition: ERF_ShocInterface.H:606
view_2d inv_qc_relvar
Definition: ERF_ShocInterface.H:658
view_1d surf_evap
Definition: ERF_ShocInterface.H:634
view_2d w_sec
Definition: ERF_ShocInterface.H:660
view_1d phis
Definition: ERF_ShocInterface.H:644
view_2d shoc_evap
Definition: ERF_ShocInterface.H:671
view_3d horiz_wind
Definition: ERF_ShocInterface.H:648
void check_flux_state_consistency(const double dt)
Definition: ERF_ShocInterface.cpp:954
typename SHF::view_1d< Real > view_1d
Definition: ERF_ShocInterface.H:47
amrex::MultiFab * m_tau13
Definition: ERF_ShocInterface.H:584
amrex::MultiFab * m_cons
Definition: ERF_ShocInterface.H:573
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:683
ekat::KokkosTypes< KokkosDefaultDevice > KT
Definition: ERF_ShocInterface.H:40
view_1d ustar
Definition: ERF_ShocInterface.H:662
scream::PhysicsFunctions< KokkosDefaultDevice > PF
Definition: ERF_ShocInterface.H:38
void finalize_impl()
Definition: ERF_ShocInterface.cpp:921
view_2d w3
Definition: ERF_ShocInterface.H:679
void initialize_impl()
Definition: ERF_ShocInterface.cpp:687
ekat::Unmanaged< typename KT::template view_2d< ScalarT > > uview_2d
Definition: ERF_ShocInterface.H:62
amrex::MultiFab * m_hfx3
Definition: ERF_ShocInterface.H:586
ekat::WorkspaceManager< Spack, KT::Device > WSM
Definition: ERF_ShocInterface.H:57
typename SHF::view_3d< const Spack > view_3d_const
Definition: ERF_ShocInterface.H:54
bool m_first_step
Definition: ERF_ShocInterface.H:596
Int m_num_tracers
Definition: ERF_ShocInterface.H:556
std::string name() const
Definition: ERF_ShocInterface.H:117
int m_step
Definition: ERF_ShocInterface.H:564
sview_2d surf_mom_flux
Definition: ERF_ShocInterface.H:633
typename SHF::view_2d< SHF::Spack > view_2d
Definition: ERF_ShocInterface.H:49
int m_lev
Definition: ERF_ShocInterface.H:561
SHOCPostprocess shoc_postprocess
Definition: ERF_ShocInterface.H:613
typename SHF::Spack Spack
Definition: ERF_ShocInterface.H:43
SHF::SHOCInputOutput input_output
Definition: ERF_ShocInterface.H:603
Buffer m_buffer
Definition: ERF_ShocInterface.H:599
void mf_to_kokkos_buffers()
Definition: ERF_ShocInterface.cpp:335
view_2d tk
Definition: ERF_ShocInterface.H:650
view_2d qv
Definition: ERF_ShocInterface.H:636
view_2d wthl_sec
Definition: ERF_ShocInterface.H:673
amrex::BoxArray m_ba
Definition: ERF_ShocInterface.H:570
view_2d tke
Definition: ERF_ShocInterface.H:652
view_1d tot_buff_view
Definition: ERF_ShocInterface.H:620
view_1d water_flux
Definition: ERF_ShocInterface.H:688
view_2d T_mid
Definition: ERF_ShocInterface.H:635
bool column_conservation_check
Definition: ERF_ShocInterface.H:627
void dealloc_buffers()
Definition: ERF_ShocInterface.cpp:261
Int hdtime
Definition: ERF_ShocInterface.H:558
void apply_turbulent_mountain_stress()
Definition: ERF_ShocInterface.cpp:934
typename SHF::view_3d< Spack > view_3d
Definition: ERF_ShocInterface.H:53
view_2d isotropy
Definition: ERF_ShocInterface.H:669
typename SHF::IntSmallPack IntSmallPack
Definition: ERF_ShocInterface.H:44
amrex::Real * m_w_subsid
Definition: ERF_ShocInterface.H:581
amrex::MultiFab * m_xvel
Definition: ERF_ShocInterface.H:576
amrex::Vector< int > m_col_offsets
Definition: ERF_ShocInterface.H:548
SHF::SHOCHistoryOutput history_output
Definition: ERF_ShocInterface.H:605
typename KT::template view_2d< const Real > sview_2d_const
Definition: ERF_ShocInterface.H:52
@ tabs
Definition: ERF_Kessler.H:24
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
Definition: ERF_ShocInterface.H:440
uview_2d< Spack > w3
Definition: ERF_ShocInterface.H:504
uview_2d< Spack > qc_copy
Definition: ERF_ShocInterface.H:491
uview_2d< Spack > dse
Definition: ERF_ShocInterface.H:489
uview_1d< Real > upwp_sfc
Definition: ERF_ShocInterface.H:458
uview_2d< Spack > vw_sec
Definition: ERF_ShocInterface.H:503
uview_1d< Real > wpthlp_sfc
Definition: ERF_ShocInterface.H:456
uview_2d< Spack > thv
Definition: ERF_ShocInterface.H:480
uview_2d< Spack > wqls_sec
Definition: ERF_ShocInterface.H:505
uview_2d< Spack > wthl_sec
Definition: ERF_ShocInterface.H:499
uview_2d< Spack > thlm
Definition: ERF_ShocInterface.H:487
uview_2d< Spack > uw_sec
Definition: ERF_ShocInterface.H:502
uview_2d< Spack > z_mid
Definition: ERF_ShocInterface.H:476
uview_1d< Real > vpwp_sfc
Definition: ERF_ShocInterface.H:459
uview_2d< Spack > rrho
Definition: ERF_ShocInterface.H:478
uview_2d< Spack > wtke_sec
Definition: ERF_ShocInterface.H:501
uview_2d< Spack > rrho_i
Definition: ERF_ShocInterface.H:479
uview_2d< Spack > unused
Definition: ERF_ShocInterface.H:475
Spack * wsm_data
Definition: ERF_ShocInterface.H:516
static constexpr int num_2d_vector_int
Definition: ERF_ShocInterface.H:449
uview_2d< Spack > wtracer_sfc
Definition: ERF_ShocInterface.H:484
static constexpr int num_2d_vector_tr
Definition: ERF_ShocInterface.H:454
static constexpr int num_1d_scalar_nlev
Definition: ERF_ShocInterface.H:446
uview_2d< Spack > inv_exner
Definition: ERF_ShocInterface.H:486
uview_2d< Spack > qwthl_sec
Definition: ERF_ShocInterface.H:498
uview_2d< Spack > qw
Definition: ERF_ShocInterface.H:488
uview_2d< Spack > brunt
Definition: ERF_ShocInterface.H:506
uview_2d< Spack > shoc_mix
Definition: ERF_ShocInterface.H:493
uview_2d< Spack > tke_copy
Definition: ERF_ShocInterface.H:490
static constexpr int num_1d_scalar_ncol
Definition: ERF_ShocInterface.H:442
uview_2d< Spack > isotropy
Definition: ERF_ShocInterface.H:494
uview_2d< Spack > zt_grid
Definition: ERF_ShocInterface.H:482
uview_2d< Spack > w_sec
Definition: ERF_ShocInterface.H:495
uview_2d< Spack > thl_sec
Definition: ERF_ShocInterface.H:496
uview_2d< Spack > wm_zt
Definition: ERF_ShocInterface.H:485
uview_2d< Spack > qw_sec
Definition: ERF_ShocInterface.H:497
uview_2d< Spack > wqw_sec
Definition: ERF_ShocInterface.H:500
uview_2d< Spack > z_int
Definition: ERF_ShocInterface.H:477
uview_2d< Spack > shoc_ql2
Definition: ERF_ShocInterface.H:492
uview_1d< Real > wprtp_sfc
Definition: ERF_ShocInterface.H:457
uview_2d< Spack > zi_grid
Definition: ERF_ShocInterface.H:483
static constexpr int num_2d_vector_mid
Definition: ERF_ShocInterface.H:448
uview_1d< Spack > pref_mid
Definition: ERF_ShocInterface.H:473
uview_2d< Spack > dz
Definition: ERF_ShocInterface.H:481
Definition: ERF_ShocInterface.H:326
view_1d_const surf_sens_flux
Definition: ERF_ShocInterface.H:391
view_1d ice_flux
Definition: ERF_ShocInterface.H:394
view_1d vapor_flux
Definition: ERF_ShocInterface.H:392
view_1d_const phis
Definition: ERF_ShocInterface.H:388
int nlev
Definition: ERF_ShocInterface.H:378
view_1d heat_flux
Definition: ERF_ShocInterface.H:395
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:384
view_2d_const rrho
Definition: ERF_ShocInterface.H:379
KOKKOS_INLINE_FUNCTION void operator()(const Kokkos::TeamPolicy< KT::ExeSpace >::member_type &team) const
Definition: ERF_ShocInterface.H:330
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:382
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:398
view_2d_const qw
Definition: ERF_ShocInterface.H:381
bool compute_mass_and_energy_fluxes
Definition: ERF_ShocInterface.H:389
view_1d water_flux
Definition: ERF_ShocInterface.H:393
view_1d_const surf_evap
Definition: ERF_ShocInterface.H:390
view_2d_const qc_copy
Definition: ERF_ShocInterface.H:381
view_2d_const z_mid
Definition: ERF_ShocInterface.H:387
view_2d_const dse
Definition: ERF_ShocInterface.H:387
view_2d tke
Definition: ERF_ShocInterface.H:380
view_2d T_mid
Definition: ERF_ShocInterface.H:386
view_2d inv_qc_relvar
Definition: ERF_ShocInterface.H:385
view_2d qc
Definition: ERF_ShocInterface.H:380
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:424
view_2d qv
Definition: ERF_ShocInterface.H:380
int ncol
Definition: ERF_ShocInterface.H:378
view_2d_const qc2
Definition: ERF_ShocInterface.H:383
view_2d_const tke_copy
Definition: ERF_ShocInterface.H:381
Definition: ERF_ShocInterface.H:129
int ncol
Definition: ERF_ShocInterface.H:223
view_3d_strided qtracers
Definition: ERF_ShocInterface.H:234
int nlev
Definition: ERF_ShocInterface.H:223
view_2d cldfrac_liq_prev
Definition: ERF_ShocInterface.H:260
view_2d rrho
Definition: ERF_ShocInterface.H:243
view_2d_const pseudo_density
Definition: ERF_ShocInterface.H:228
sview_2d_const surf_mom_flux
Definition: ERF_ShocInterface.H:233
view_2d_const qc
Definition: ERF_ShocInterface.H:236
view_2d thlm
Definition: ERF_ShocInterface.H:256
view_1d_const phis
Definition: ERF_ShocInterface.H:230
view_1d_const surf_evap
Definition: ERF_ShocInterface.H:232
view_2d zt_grid
Definition: ERF_ShocInterface.H:247
view_2d cldfrac_liq
Definition: ERF_ShocInterface.H:259
view_2d inv_exner
Definition: ERF_ShocInterface.H:255
view_1d_const surf_sens_flux
Definition: ERF_ShocInterface.H:231
view_2d wm_zt
Definition: ERF_ShocInterface.H:254
KOKKOS_INLINE_FUNCTION void operator()(const Kokkos::TeamPolicy< KT::ExeSpace >::member_type &team) const
Definition: ERF_ShocInterface.H:133
view_2d qc_copy
Definition: ERF_ShocInterface.H:237
view_2d wtracer_sfc
Definition: ERF_ShocInterface.H:253
view_2d z_mid
Definition: ERF_ShocInterface.H:238
view_2d rrho_i
Definition: ERF_ShocInterface.H:244
view_2d_const omega
Definition: ERF_ShocInterface.H:229
view_2d_const T_mid
Definition: ERF_ShocInterface.H:225
view_2d cloud_frac
Definition: ERF_ShocInterface.H:258
Real z_surf
Definition: ERF_ShocInterface.H:224
view_2d_const p_int
Definition: ERF_ShocInterface.H:227
view_2d tke
Definition: ERF_ShocInterface.H:241
view_1d wprtp_sfc
Definition: ERF_ShocInterface.H:250
view_2d qw
Definition: ERF_ShocInterface.H:257
view_1d upwp_sfc
Definition: ERF_ShocInterface.H:251
view_1d wpthlp_sfc
Definition: ERF_ShocInterface.H:249
view_2d z_int
Definition: ERF_ShocInterface.H:239
view_2d_const p_mid
Definition: ERF_ShocInterface.H:226
view_1d vpwp_sfc
Definition: ERF_ShocInterface.H:252
view_2d dz
Definition: ERF_ShocInterface.H:246
view_2d shoc_s
Definition: ERF_ShocInterface.H:240
view_2d tke_copy
Definition: ERF_ShocInterface.H:242
view_2d zi_grid
Definition: ERF_ShocInterface.H:248
view_2d qv
Definition: ERF_ShocInterface.H:235
view_2d thv
Definition: ERF_ShocInterface.H:245
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:263
Definition: ERF_DataStruct.H:123