ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SlowRhsPre.cpp File Reference
#include "AMReX_MultiFab.H"
#include "AMReX_iMultiFab.H"
#include "AMReX_ArrayLim.H"
#include "AMReX_BCRec.H"
#include "AMReX_GpuContainers.H"
#include "AMReX_GpuPrint.H"
#include "ERF_TI_slow_headers.H"
#include "ERF_EOS.H"
#include "ERF_Utils.H"
#include "ERF_Diffusion.H"
#include "ERF_EBAdvection.H"
#include "ERF_EB.H"
#include "ERF_SurfaceLayer.H"
Include dependency graph for ERF_SlowRhsPre.cpp:

Functions

void erf_slow_rhs_pre (int level, int finest_level, int nrk, Real dt, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old, Vector< MultiFab > &S_data, const MultiFab &S_prim, const MultiFab &qt, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &zvel, std::unique_ptr< MultiFab > &z_t_mf, const MultiFab &cc_src, const MultiFab &xmom_src, const MultiFab &ymom_src, const MultiFab &zmom_src, const MultiFab &buoyancy, const MultiFab *zmom_crse_rhs, Vector< std::unique_ptr< MultiFab >> &Tau_lev, Vector< std::unique_ptr< MultiFab >> &Tau_corr_lev, Vector< std::unique_ptr< MultiFab >> &Tau_EB, MultiFab *SmnSmn, MultiFab *eddyDiffs, MultiFab *Hfx1, MultiFab *Hfx2, MultiFab *Hfx3, MultiFab *Q1fx1, MultiFab *Q1fx2, MultiFab *Q1fx3, MultiFab *Q2fx3, MultiFab *Diss, MultiFab *Hfx3_EB, const Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, const Gpu::DeviceVector< BCRec > &domain_bcs_type_d, const Vector< BCRec > &domain_bcs_type_h, const MultiFab &z_phys_nd, const MultiFab &z_phys_cc, const MultiFab &ax, const MultiFab &ay, const MultiFab &az, const MultiFab &detJ, Gpu::DeviceVector< Real > &stretched_dz_d, Vector< MultiFab > &gradp, Vector< std::unique_ptr< MultiFab >> &mapfac, const eb_ &ebfact, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine)
 

Function Documentation

◆ erf_slow_rhs_pre()

void erf_slow_rhs_pre ( int  level,
int  finest_level,
int  nrk,
Real  dt,
Vector< MultiFab > &  S_rhs,
Vector< MultiFab > &  S_old,
Vector< MultiFab > &  S_data,
const MultiFab &  S_prim,
const MultiFab &  qt,
MultiFab &  avg_xmom,
MultiFab &  avg_ymom,
MultiFab &  avg_zmom,
const MultiFab &  xvel,
const MultiFab &  yvel,
const MultiFab &  zvel,
std::unique_ptr< MultiFab > &  z_t_mf,
const MultiFab &  cc_src,
const MultiFab &  xmom_src,
const MultiFab &  ymom_src,
const MultiFab &  zmom_src,
const MultiFab &  buoyancy,
const MultiFab *  zmom_crse_rhs,
Vector< std::unique_ptr< MultiFab >> &  Tau_lev,
Vector< std::unique_ptr< MultiFab >> &  Tau_corr_lev,
Vector< std::unique_ptr< MultiFab >> &  Tau_EB,
MultiFab *  SmnSmn,
MultiFab *  eddyDiffs,
MultiFab *  Hfx1,
MultiFab *  Hfx2,
MultiFab *  Hfx3,
MultiFab *  Q1fx1,
MultiFab *  Q1fx2,
MultiFab *  Q1fx3,
MultiFab *  Q2fx3,
MultiFab *  Diss,
MultiFab *  Hfx3_EB,
const Geometry  geom,
const SolverChoice solverChoice,
std::unique_ptr< SurfaceLayer > &  SurfLayer,
const Gpu::DeviceVector< BCRec > &  domain_bcs_type_d,
const Vector< BCRec > &  domain_bcs_type_h,
const MultiFab &  z_phys_nd,
const MultiFab &  z_phys_cc,
const MultiFab &  ax,
const MultiFab &  ay,
const MultiFab &  az,
const MultiFab &  detJ,
Gpu::DeviceVector< Real > &  stretched_dz_d,
Vector< MultiFab > &  gradp,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
const eb_ ebfact,
YAFluxRegister *  fr_as_crse,
YAFluxRegister *  fr_as_fine 
)

Function for computing the slow RHS for the evolution equations for the density, potential temperature and momentum.

Parameters
[in]levellevel of resolution
[in]finest_levelfinest level of resolution
[in]nrkwhich RK stage
[in]dtslow time step
[out]S_rhsRHS computed here
[in]S_oldold-time solution – used only for anelastic
[in]S_datacurrent solution
[in]S_primprimitive variables (i.e. conserved variables divided by density)
[in,out]avg_xmom
[in,out]avg_ymom
[in,out]avg_zmom
[in]xvelx-component of velocity
[in]yvely-component of velocity
[in]zvelz-component of velocity
[in]z_t_mf rate of change of grid height – only relevant for moving terrain
[in]cc_srcsource terms for conserved variables
[in]xmom_srcsource terms for x-momentum
[in]ymom_srcsource terms for y-momentum
[in]zmom_srcsource terms for z-momentum
[in]buoyancybuoyancy source term for z-momentum
[in]zmom_crse_rhsupdate term from coarser level for z-momentum; non-zero on c/f boundary only
[in]Tau_levcomponents of stress tensor
[in]SmnSmnstrain rate magnitude
[in]eddyDiffsdiffusion coefficients for LES turbulence models
[in]Hfx3heat flux in z-dir
[in]Dissdissipation of turbulent kinetic energy
[in]geomContainer for geometric information
[in]solverChoiceContainer for solver parameters
[in]SurfLayerPointer to SurfaceLayer class for Monin-Obukhov Similarity Theory boundary condition
[in]domain_bcs_type_ddevice vector for domain boundary conditions
[in]domain_bcs_type_hhost vector for domain boundary conditions
[in]z_phys_ndheight coordinate at nodes
[in]axarea fractions on x-faces
[in]ayarea fractions on y-faces
[in]azarea fractions on z-faces
[in]detJJacobian of the metric transformation (= 1 if use_terrain_fitted_coords is false)
[in]gradppressure gradient
[in]mapfacmap factors
[in]ebfactEB factories for cell- and face-centered variables
[in,out]fr_as_crseYAFluxRegister at level l at level l / l+1 interface
[in,out]fr_as_fineYAFluxRegister at level l at level l-1 / l interface
114 {
115  BL_PROFILE_REGION("erf_slow_rhs_pre()");
116 
117  const BCRec* bc_ptr_d = domain_bcs_type_d.data();
118  const BCRec* bc_ptr_h = domain_bcs_type_h.data();
119 
120  DiffChoice dc = solverChoice.diffChoice;
121  TurbChoice tc = solverChoice.turbChoice[level];
122 
123  const MultiFab* t_mean_mf = nullptr;
124  if (SurfLayer) { t_mean_mf = SurfLayer->get_mac_avg(level,2); }
125 
126  const Box& domain = geom.Domain();
127  int klo = domain.smallEnd(2);
128  int khi = domain.bigEnd(2);
129 
130  const AdvType l_horiz_adv_type = solverChoice.advChoice.dycore_horiz_adv_type;
131  const AdvType l_vert_adv_type = solverChoice.advChoice.dycore_vert_adv_type;
132  const Real l_horiz_upw_frac = solverChoice.advChoice.dycore_horiz_upw_frac;
133  const Real l_vert_upw_frac = solverChoice.advChoice.dycore_vert_upw_frac;
134  const bool l_use_stretched_dz = (solverChoice.mesh_type == MeshType::StretchedDz);
135  const bool l_use_terrain_fitted_coords = (solverChoice.mesh_type == MeshType::VariableDz);
136  const bool l_moving_terrain = (solverChoice.terrain_type == TerrainType::MovingFittedMesh);
137  if (l_moving_terrain) AMREX_ALWAYS_ASSERT (l_use_stretched_dz || l_use_terrain_fitted_coords);
138 
139  const bool l_use_diff = ( (dc.molec_diff_type != MolecDiffType::None) ||
140  (tc.les_type != LESType::None) ||
141  (tc.rans_type != RANSType::None) ||
142  (tc.pbl_type != PBLType::None) );
143  const bool l_use_turb = tc.use_kturb;
144  const bool l_need_SmnSmn = tc.use_keqn;
145 
146  const Real l_vert_implicit_fac = (solverChoice.implicit_thermal_diffusion) ?
147  solverChoice.vert_implicit_fac[nrk] : 0.0;
148 
149  const bool l_use_moisture = (solverChoice.moisture_type != MoistureType::None);
150  const bool l_use_SurfLayer = (SurfLayer != nullptr);
151  const bool l_rotate = (solverChoice.use_rotate_surface_flux);
152 
153  const bool l_anelastic = (solverChoice.anelastic[level] == 1);
154  const bool l_fixed_rho = (solverChoice.fixed_density[level] == 1);
155 
156  const bool l_reflux = ( (solverChoice.coupling_type == CouplingType::TwoWay) && (finest_level > 0) &&
157  ( (l_anelastic && nrk == 1) || (!l_anelastic && nrk == 2) ) );
158 
159  const bool l_use_eb = (solverChoice.terrain_type == TerrainType::EB);
160 
161  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
162  const Real* dx = geom.CellSize();
163 
164  // *****************************************************************************
165  // Combine external forcing terms
166  // *****************************************************************************
167  const Array<Real,AMREX_SPACEDIM> grav{zero, zero, -solverChoice.gravity};
168  const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
169 
170  // **************************************************************************************
171  // If doing advection with EB we need the extra values for tangential interpolation
172  // **************************************************************************************
173  if (l_use_eb) {
174  S_data[IntVars::xmom].FillBoundary(geom.periodicity());
175  S_data[IntVars::ymom].FillBoundary(geom.periodicity());
176  S_data[IntVars::zmom].FillBoundary(geom.periodicity());
177  }
178 
179  // *****************************************************************************
180  // Pre-computed quantities
181  // *****************************************************************************
182  int nvars = S_data[IntVars::cons].nComp();
183  const BoxArray& ba = S_data[IntVars::cons].boxArray();
184  const DistributionMapping& dm = S_data[IntVars::cons].DistributionMap();
185 
186  int nGhost = (l_use_eb) ? 2 : 1;
187  MultiFab Omega(convert(ba,IntVect(0,0,1)), dm, 1, nGhost);
188 
189  std::unique_ptr<MultiFab> expr;
190  std::unique_ptr<MultiFab> dflux_x;
191  std::unique_ptr<MultiFab> dflux_y;
192  std::unique_ptr<MultiFab> dflux_z;
193 
194  if (l_use_diff) {
195 #ifdef ERF_USE_SHOC
196  if (solverChoice.use_shoc) {
197  // Populate vertical component of eddyDiffs
198  shoc_lev->set_eddy_diffs();
199  }
200 #endif
201 
202  erf_make_tau_terms(level,nrk,domain_bcs_type_h,z_phys_nd,
203  S_data,xvel,yvel,zvel,
204  Tau_lev,Tau_corr_lev,
205  SmnSmn,eddyDiffs,geom,solverChoice,SurfLayer,
206  stretched_dz_d, detJ,mapfac, ax, ay, az, ebfact);
207 
208  IntVect ng(0,0,1);
209  dflux_x = std::make_unique<MultiFab>(convert(ba,IntVect(1,0,0)), dm, nvars, ng);
210  dflux_y = std::make_unique<MultiFab>(convert(ba,IntVect(0,1,0)), dm, nvars, ng);
211  dflux_z = std::make_unique<MultiFab>(convert(ba,IntVect(0,0,1)), dm, nvars, 0);
212 
213 #ifdef ERF_USE_SHOC
214  if (solverChoice.use_shoc) {
215  // Zero out the surface stresses of tau13/tau23
216  shoc_lev->set_diff_stresses();
217  } else if (l_use_SurfLayer) {
218  // Set surface shear stresses, update heat and moisture fluxes
219  // (fluxes will be later applied in the diffusion source update)
220  Vector<const MultiFab*> mfs = {&S_data[IntVars::cons], &xvel, &yvel, &zvel};
221  SurfLayer->impose_SurfaceLayer_bcs(level, mfs, Tau_lev,
222  Hfx1, Hfx2, Hfx3,
223  Q1fx1, Q1fx2, Q1fx3,
224  &z_phys_nd);
225  }
226 #else
227  // This is computed pre step in Advance if we use SHOC
228  if (l_use_SurfLayer) {
229  if (!l_use_eb) {
230  // Set surface shear stresses, update heat and moisture fluxes
231  // (fluxes will be later applied in the diffusion source update)
232  Vector<const MultiFab*> mfs = {&S_data[IntVars::cons], &xvel, &yvel, &zvel};
233  SurfLayer->impose_SurfaceLayer_bcs(level, mfs, Tau_lev,
234  Hfx1, Hfx2, Hfx3,
235  Q1fx1, Q1fx2, Q1fx3,
236  &z_phys_nd);
237 
238  //if (l_vert_implicit_fac > 0 && solverChoice.implicit_momentum_diffusion) {
239  // copy_surface_tau_for_implicit(Tau_lev, Tau_corr_lev);
240  //}
241  } else {
242  Vector<const MultiFab*> mfs = {&S_data[IntVars::cons], &xvel, &yvel, &zvel};
243  SurfLayer->impose_SurfaceLayer_bcs_EB(level, mfs, Tau_EB,
244  Hfx1, Hfx2, Hfx3_EB,
245  Q1fx1, Q1fx2, Q1fx3,
246  ebfact);
247  }
248  }
249 #endif
250  } // l_use_diff
251 
252  // This is just cautionary to deal with grid boundaries that aren't domain boundaries
253  S_rhs[IntVars::zmom].setVal(0);
254 
255  // *****************************************************************************
256  // Define updates and fluxes in the current RK stage
257  // *****************************************************************************
258  // Cell-centered masks for EB (used for flux interpolation)
259  bool already_on_centroids = false;
260  Vector<iMultiFab> physbnd_mask;
261  physbnd_mask.resize(IntVars::NumTypes);
262  if (l_use_eb) {
263  physbnd_mask[IntVars::cons].define(S_data[IntVars::cons].boxArray(), S_data[IntVars::cons].DistributionMap(), 1, 1);
264  physbnd_mask[IntVars::cons].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
265  // physbnd_mask[IntVars::cons].FillBoundary(geom.periodicity());
266  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
267  physbnd_mask[1+dir].define(S_data[1+dir].boxArray(), S_data[1+dir].DistributionMap(), 1, 1);
268  physbnd_mask[1+dir].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
269  // physbnd_mask[1+dir].FillBoundary(geom.periodicity());
270  }
271  }
272 
273 #ifdef _OPENMP
274 #pragma omp parallel if (Gpu::notInLaunchRegion())
275 #endif
276  {
277  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
278  {
279  Box bx = mfi.tilebox();
280  Box tbx = mfi.nodaltilebox(0);
281  Box tby = mfi.nodaltilebox(1);
282  Box tbz = mfi.nodaltilebox(2);
283 
284  // Boxes for momentum fluxes
285  Vector<Box> tbx_grown(AMREX_SPACEDIM);
286  Vector<Box> tby_grown(AMREX_SPACEDIM);
287  Vector<Box> tbz_grown(AMREX_SPACEDIM);
288  if (l_use_eb) {
289  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
290  tbx_grown[dir] = tbx;
291  tby_grown[dir] = tby;
292  tbz_grown[dir] = tbz;
293  IntVect iv(1, 1, 1);
294  iv[dir] = 0;
295  tbx_grown[dir] = (tbx_grown[dir].growHi(dir,1)).grow(iv);
296  tby_grown[dir] = (tby_grown[dir].growHi(dir,1)).grow(iv);
297  tbz_grown[dir] = (tbz_grown[dir].growHi(dir,1)).grow(iv);
298  }
299  }
300 
301  // We don't compute a source term for z-momentum on the bottom or top domain boundary
302  if (tbz.smallEnd(2) == domain.smallEnd(2)) {
303  tbz.growLo(2,-1);
304  }
305  if (tbz.bigEnd(2) == domain.bigEnd(2)+1) {
306  tbz.growHi(2,-1);
307  }
308 
309  const Array4<const Real> & cell_data = S_data[IntVars::cons].array(mfi);
310  const Array4<const Real> & cell_prim = S_prim.array(mfi);
311  const Array4<Real> & cell_rhs = S_rhs[IntVars::cons].array(mfi);
312 
313  const Array4<const Real> & cell_old = S_old[IntVars::cons].array(mfi);
314 
315  const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
316  const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
317  const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
318  const Array4<Real const>& buoyancy_arr = buoyancy.const_array(mfi);
319 
320  const Array4<Real const>& gpx_arr = gradp[GpVars::gpx].const_array(mfi);
321  const Array4<Real const>& gpy_arr = gradp[GpVars::gpy].const_array(mfi);
322  const Array4<Real const>& gpz_arr = gradp[GpVars::gpz].const_array(mfi);
323 
324  const Array4<Real const>& qt_arr = qt.const_array(mfi);
325 
326  const Array4<Real>& rho_u_old = S_old[IntVars::xmom].array(mfi);
327  const Array4<Real>& rho_v_old = S_old[IntVars::ymom].array(mfi);
328 
329  if (l_anelastic) {
330  // When anelastic we must reset these to 0 each RK step
331  avg_xmom[mfi].template setVal<RunOn::Device>(0,tbx);
332  avg_ymom[mfi].template setVal<RunOn::Device>(0,tby);
333  avg_zmom[mfi].template setVal<RunOn::Device>(0,tbz);
334  }
335 
336  Array4<Real> avg_xmom_arr = avg_xmom.array(mfi);
337  Array4<Real> avg_ymom_arr = avg_ymom.array(mfi);
338  Array4<Real> avg_zmom_arr = avg_zmom.array(mfi);
339 
340  const Array4<const Real> & u = xvel.array(mfi);
341  const Array4<const Real> & v = yvel.array(mfi);
342  const Array4<const Real> & w = zvel.array(mfi);
343 
344  const Array4<const Real>& rho_u = S_data[IntVars::xmom].array(mfi);
345  const Array4<const Real>& rho_v = S_data[IntVars::ymom].array(mfi);
346  const Array4<const Real>& rho_w = S_data[IntVars::zmom].array(mfi);
347 
348  // Map factors
349  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
350  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
351  const Array4<const Real>& mf_vx = mapfac[MapFacType::v_x]->const_array(mfi);
352  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
353  const Array4<const Real>& mf_uy = mapfac[MapFacType::u_y]->const_array(mfi);
354  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
355 
356  const Array4< Real>& omega_arr = Omega.array(mfi);
357 
358  Array4<const Real> z_t;
359  if (z_t_mf) {
360  z_t = z_t_mf->array(mfi);
361  } else {
362  z_t = Array4<const Real>{};
363  }
364 
365  const Array4<Real>& rho_u_rhs = S_rhs[IntVars::xmom].array(mfi);
366  const Array4<Real>& rho_v_rhs = S_rhs[IntVars::ymom].array(mfi);
367  const Array4<Real>& rho_w_rhs = S_rhs[IntVars::zmom].array(mfi);
368 
369  const Array4<Real const>& mu_turb = l_use_turb ? eddyDiffs->const_array(mfi) : Array4<const Real>{};
370 
371  // Terrain metrics
372  const Array4<const Real>& z_nd = z_phys_nd.const_array(mfi);
373  const Array4<const Real>& z_cc = z_phys_cc.const_array(mfi);
374 
375  // *****************************************************************************
376  // Define flux arrays for use in advection
377  // *****************************************************************************
378  std::array<FArrayBox,AMREX_SPACEDIM> flux;
379  std::array<FArrayBox,AMREX_SPACEDIM> flux_u;
380  std::array<FArrayBox,AMREX_SPACEDIM> flux_v;
381  std::array<FArrayBox,AMREX_SPACEDIM> flux_w;
382 
383  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
384  if (!l_use_eb) {
385  flux[dir].resize(surroundingNodes(bx,dir),2,The_Async_Arena());
386  } else {
387  flux[dir].resize(surroundingNodes(bx,dir).grow(1),2,The_Async_Arena());
388  }
389  flux[dir].setVal<RunOn::Device>(0);
390  }
391  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
392  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
393 
394  // Define flux arrays for momentum variables (used only for EB now)
395  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_u_arr{};
396  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_v_arr{};
397  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_w_arr{};
398 
399  if (l_use_eb) {
400  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
401  flux_u[dir].resize(tbx_grown[dir],1,The_Async_Arena());
402  flux_v[dir].resize(tby_grown[dir],1,The_Async_Arena());
403  flux_w[dir].resize(tbz_grown[dir],1,The_Async_Arena());
404  flux_u[dir].setVal<RunOn::Device>(0);
405  flux_v[dir].setVal<RunOn::Device>(0);
406  flux_w[dir].setVal<RunOn::Device>(0);
407  flx_u_arr[dir] = flux_u[dir].array();
408  flx_v_arr[dir] = flux_v[dir].array();
409  flx_w_arr[dir] = flux_w[dir].array();
410  }
411  }
412 
413  // *****************************************************************************
414  // Contravariant flux field
415  // *****************************************************************************
416  {
417  BL_PROFILE("slow_rhs_making_omega");
418  IntVect nGrowVect = (l_use_eb)
419  ? IntVect(AMREX_D_DECL(2, 2, 2)) : IntVect(AMREX_D_DECL(1, 1, 1));
420  Box gbxo = surroundingNodes(bx,2); gbxo.grow(nGrowVect);
421  //
422  // Now create Omega with momentum (not velocity) with z_t subtracted if moving terrain
423  // ONLY if not doing anelastic + terrain -- in that case Omega will be defined coming
424  // out of the projection
425  //
426  if (!l_use_terrain_fitted_coords) {
427  ParallelFor(gbxo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
428  omega_arr(i,j,k) = rho_w(i,j,k);
429  });
430 
431  } else {
432 
433  Box gbxo_lo = gbxo; gbxo_lo.setBig(2,domain.smallEnd(2));
434  int lo_z_face = domain.smallEnd(2);
435  if (gbxo_lo.smallEnd(2) <= lo_z_face) {
436  ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
437  omega_arr(i,j,k) = zero;
438  });
439  }
440  Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
441  int hi_z_face = domain.bigEnd(2)+1;
442  if (gbxo_hi.bigEnd(2) >= hi_z_face) {
443  ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
444  omega_arr(i,j,k) = rho_w(i,j,k);
445  });
446  }
447 
448  if (z_t) { // Note we never do anelastic with moving terrain
449  Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
450  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
451  // We define rho on the z-face the same way as in MomentumToVelocity/VelocityToMomentum
452  Real rho_at_face = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp));
453  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
454  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv) -
455  rho_at_face * z_t(i,j,k);
456  });
457  } else {
458  Box gbxo_mid = gbxo;
459  if (gbxo_mid.smallEnd(2) <= domain.smallEnd(2)) {
460  gbxo_mid.setSmall(2,1);
461  }
462  if (gbxo_mid.bigEnd(2) >= domain.bigEnd(2)+1) {
463  gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
464  }
465  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
466  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
467  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv);
468  });
469  }
470  }
471  } // end profile
472 
473  // *****************************************************************************
474  // Diffusive terms (pre-computed above)
475  // *****************************************************************************
476  // No terrain diffusion
477  Array4<Real> tau11,tau22,tau33;
478  Array4<Real> tau12,tau13,tau23;
479  if (Tau_lev[TauType::tau11]) {
480  tau11 = Tau_lev[TauType::tau11]->array(mfi); tau22 = Tau_lev[TauType::tau22]->array(mfi);
481  tau33 = Tau_lev[TauType::tau33]->array(mfi); tau12 = Tau_lev[TauType::tau12]->array(mfi);
482  tau13 = Tau_lev[TauType::tau13]->array(mfi); tau23 = Tau_lev[TauType::tau23]->array(mfi);
483  } else {
484  tau11 = Array4<Real>{}; tau22 = Array4<Real>{}; tau33 = Array4<Real>{};
485  tau12 = Array4<Real>{}; tau13 = Array4<Real>{}; tau23 = Array4<Real>{};
486  }
487  // Terrain diffusion
488  Array4<Real> tau21,tau31,tau32;
489  if (Tau_lev[TauType::tau21]) {
490  tau21 = Tau_lev[TauType::tau21]->array(mfi);
491  tau31 = Tau_lev[TauType::tau31]->array(mfi);
492  tau32 = Tau_lev[TauType::tau32]->array(mfi);
493  } else {
494  tau21 = Array4<Real>{}; tau31 = Array4<Real>{}; tau32 = Array4<Real>{};
495  }
496 
497  Array4<Real> tau_eb13{}, tau_eb23{};
498  if (l_use_eb && Tau_EB[EBTauType::tau_eb13] && Tau_EB[EBTauType::tau_eb23]) {
499  tau_eb13 = Tau_EB[EBTauType::tau_eb13]->array(mfi);
500  tau_eb23 = Tau_EB[EBTauType::tau_eb23]->array(mfi);
501  }
502 
503  // Strain magnitude
504  Array4<Real> SmnSmn_a;
505  if (l_need_SmnSmn) {
506  SmnSmn_a = SmnSmn->array(mfi);
507  } else {
508  SmnSmn_a = Array4<Real>{};
509  }
510 
511  // *****************************************************************************
512  // Define updates in the RHS of continuity and potential temperature equations
513  // *****************************************************************************
514  bool l_eb_terrain_cc = false; // EB terrain on cell-centered grid
515  Array4<const int> mask_arr{};
516  Array4<const EBCellFlag> cfg_arr{};
517  Array4<const Real> ax_arr{};
518  Array4<const Real> ay_arr{};
519  Array4<const Real> az_arr{};
520  Array4<const Real> fcx_arr{};
521  Array4<const Real> fcy_arr{};
522  Array4<const Real> fcz_arr{};
523  Array4<const Real> detJ_arr{};
524  Array4<const Real> barea_arr{};
525  Array4<const Real> bcent_arr{};
526 
527  if (l_use_eb)
528  {
529  EBCellFlagFab const& cfg = (ebfact.get_const_factory())->getMultiEBCellFlagFab()[mfi];
530  cfg_arr = cfg.const_array();
531  if (cfg.getType(bx) == FabType::singlevalued) {
532  l_eb_terrain_cc = true;
533  ax_arr = (ebfact.get_const_factory())->getAreaFrac()[0]->const_array(mfi);
534  ay_arr = (ebfact.get_const_factory())->getAreaFrac()[1]->const_array(mfi);
535  az_arr = (ebfact.get_const_factory())->getAreaFrac()[2]->const_array(mfi);
536  fcx_arr = (ebfact.get_const_factory())->getFaceCent()[0]->const_array(mfi);
537  fcy_arr = (ebfact.get_const_factory())->getFaceCent()[1]->const_array(mfi);
538  fcz_arr = (ebfact.get_const_factory())->getFaceCent()[2]->const_array(mfi);
539  detJ_arr = (ebfact.get_const_factory())->getVolFrac().const_array(mfi);
540  mask_arr = physbnd_mask[IntVars::cons].const_array(mfi);
541  barea_arr = (ebfact.get_const_factory())->getBndryArea().const_array(mfi);
542  bcent_arr = (ebfact.get_const_factory())->getBndryCent().const_array(mfi);
543  } else {
544  ax_arr = ax.const_array(mfi);
545  ay_arr = ay.const_array(mfi);
546  az_arr = az.const_array(mfi);
547  detJ_arr = detJ.const_array(mfi);
548  }
549  } else {
550  ax_arr = ax.const_array(mfi);
551  ay_arr = ay.const_array(mfi);
552  az_arr = az.const_array(mfi);
553  detJ_arr = detJ.const_array(mfi);
554  }
555 
556  int icomp = RhoTheta_comp; int ncomp = 1;
557  if (!l_eb_terrain_cc){
558  AdvectionSrcForRho( bx, cell_rhs,
559  rho_u, rho_v, omega_arr, // these are being used to build the fluxes
560  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr, // these are being defined from the fluxes
561  ax_arr, ay_arr, az_arr, detJ_arr,
562  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
563  flx_arr, l_fixed_rho);
564  AdvectionSrcForScalars(bx, icomp, ncomp,
565  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
566  cell_prim, cell_rhs,
567  detJ_arr, dxInv, mf_mx, mf_my,
568  l_horiz_adv_type, l_vert_adv_type,
569  l_horiz_upw_frac, l_vert_upw_frac,
570  flx_arr, domain, bc_ptr_h);
571  } else {
572  EBAdvectionSrcForRho(bx, cell_rhs,
573  rho_u, rho_v, omega_arr,
574  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
575  mask_arr, cfg_arr,
576  ax_arr, ay_arr, az_arr,
577  fcx_arr, fcy_arr, fcz_arr, detJ_arr,
578  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
579  flx_arr, l_fixed_rho,
580  already_on_centroids);
581  EBAdvectionSrcForScalars(bx, icomp, ncomp,
582  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
583  cell_prim, cell_rhs,
584  mask_arr, cfg_arr, ax_arr, ay_arr, az_arr,
585  fcx_arr, fcy_arr, fcz_arr,
586  detJ_arr, dxInv, mf_mx, mf_my,
587  l_horiz_adv_type, l_vert_adv_type,
588  l_horiz_upw_frac, l_vert_upw_frac,
589  flx_arr, domain, bc_ptr_h,
590  already_on_centroids);
591  }
592 
593  if (l_use_diff) {
594  Array4<Real> diffflux_x = dflux_x->array(mfi);
595  Array4<Real> diffflux_y = dflux_y->array(mfi);
596  Array4<Real> diffflux_z = dflux_z->array(mfi);
597 
598  Array4<Real> hfx_x = Hfx1->array(mfi);
599  Array4<Real> hfx_y = Hfx2->array(mfi);
600  Array4<Real> hfx_z = Hfx3->array(mfi);
601  Array4<Real> hfx_EB{};
602  if (l_use_eb) {
603  hfx_EB = Hfx3_EB->array(mfi);
604  }
605 
606  Array4<Real> q1fx_x = (Q1fx1) ? Q1fx1->array(mfi) : Array4<Real>{};
607  Array4<Real> q1fx_y = (Q1fx2) ? Q1fx2->array(mfi) : Array4<Real>{};
608  Array4<Real> q1fx_z = (Q1fx3) ? Q1fx3->array(mfi) : Array4<Real>{};
609 
610  Array4<Real> q2fx_z = (Q2fx3) ? Q2fx3->array(mfi) : Array4<Real>{};
611  Array4<Real> diss = Diss->array(mfi);
612 
613  const Array4<const Real> tm_arr = t_mean_mf ? t_mean_mf->const_array(mfi) : Array4<const Real>{};
614 
615  // NOTE: No diffusion for continuity, so n starts at one
616  int n_start = RhoTheta_comp;
617  int n_comp = 1;
618 
619  // For l_vert_implicit_fac > 0, we scale the rho*theta contribution
620  // by (1 - implicit_fac) and add in the implicit contribution with
621  // ERF_Implicit.H
622  if (l_use_stretched_dz) {
623  DiffusionSrcForState_S(bx, domain, n_start, n_comp, u, v,
624  cell_data, cell_prim, cell_rhs,
625  diffflux_x, diffflux_y, diffflux_z,
626  stretched_dz_d, dxInv, SmnSmn_a,
627  mf_mx, mf_ux, mf_vx,
628  mf_my, mf_uy, mf_vy,
629  hfx_z, q1fx_z, q2fx_z, diss,
630  mu_turb, solverChoice, level,
631  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
632  } else if (l_use_terrain_fitted_coords) {
633  DiffusionSrcForState_T(bx, domain, n_start, n_comp, l_rotate, u, v,
634  cell_data, cell_prim, cell_rhs,
635  diffflux_x, diffflux_y, diffflux_z,
636  z_nd, z_cc, ax_arr, ay_arr, az_arr, detJ_arr,
637  dxInv, SmnSmn_a,
638  mf_mx, mf_ux, mf_vx,
639  mf_my, mf_uy, mf_vy,
640  hfx_x, hfx_y, hfx_z, q1fx_x, q1fx_y, q1fx_z, q2fx_z, diss,
641  mu_turb, solverChoice, level,
642  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
643  } else if (l_use_eb) {
644  DiffusionSrcForState_EB(bx, domain, n_start, n_comp, u, v,
645  cell_data, cell_prim, cell_rhs,
646  diffflux_x, diffflux_y, diffflux_z,
647  cfg_arr, ax_arr, ay_arr, az_arr, detJ_arr,
648  barea_arr, bcent_arr,
649  dx, dxInv,
650  hfx_z, q1fx_z, q2fx_z, hfx_EB,
651  mu_turb, solverChoice, level,
652  bc_ptr_d, l_use_SurfLayer);
653  } else {
654  DiffusionSrcForState_N(bx, domain, n_start, n_comp, u, v,
655  cell_data, cell_prim, cell_rhs,
656  diffflux_x, diffflux_y, diffflux_z,
657  dxInv, SmnSmn_a,
658  mf_mx, mf_ux, mf_vx,
659  mf_my, mf_uy, mf_vy,
660  hfx_z, q1fx_z, q2fx_z, diss,
661  mu_turb, solverChoice, level,
662  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
663  }
664  }
665 
666  const Array4<Real const>& source_arr = cc_src.const_array(mfi);
667  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
668  {
669  cell_rhs(i,j,k,Rho_comp) += source_arr(i,j,k,Rho_comp);
670  cell_rhs(i,j,k,RhoTheta_comp) += source_arr(i,j,k,RhoTheta_comp);
671  });
672 
673  // If anelastic and in second RK stage, take average of old-time and new-time source
674  if ( l_anelastic && (nrk == 1) )
675  {
676  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
677  {
678  cell_rhs(i,j,k, Rho_comp) *= myhalf;
679  cell_rhs(i,j,k,RhoTheta_comp) *= myhalf;
680 
681  cell_rhs(i,j,k, Rho_comp) += myhalf / dt * (cell_data(i,j,k, Rho_comp) - cell_old(i,j,k, Rho_comp));
682  cell_rhs(i,j,k,RhoTheta_comp) += myhalf / dt * (cell_data(i,j,k,RhoTheta_comp) - cell_old(i,j,k,RhoTheta_comp));
683  });
684  }
685 
686  // *****************************************************************************
687  // Define updates in the RHS of {x, y, z}-momentum equations
688  // *****************************************************************************
689  int lo_z_face = domain.smallEnd(2);
690  int hi_z_face = domain.bigEnd(2)+1;
691 
692  AdvectionSrcForMom(mfi, bx, tbx, tby, tbz, tbx_grown, tby_grown, tbz_grown,
693  rho_u_rhs, rho_v_rhs, rho_w_rhs,
694  cell_data, u, v, w,
695  rho_u, rho_v, omega_arr,
696  z_nd, ax_arr, ay_arr, az_arr,
697  detJ_arr, stretched_dz_d,
698  dxInv, mf_mx, mf_ux, mf_vx, mf_my, mf_uy, mf_vy,
699  l_horiz_adv_type, l_vert_adv_type,
700  l_horiz_upw_frac, l_vert_upw_frac,
701  solverChoice.mesh_type, solverChoice.terrain_type,
702  ebfact, flx_u_arr, flx_v_arr, flx_w_arr,
703  physbnd_mask, already_on_centroids,
704  lo_z_face, hi_z_face, domain, bc_ptr_h);
705 
706  if (l_use_diff) {
707  // Note: tau** were calculated with calls to
708  // ComputeStress[Cons|Var]Visc_[N|S|T] in which ConsVisc ("constant
709  // viscosity") means that there is no contribution from a
710  // turbulence model. However, whether this field truly is constant
711  // depends on whether MolecDiffType is Constant or ConstantAlpha.
712  if (!l_use_eb) {
713  DiffusionSrcForMom(tbx, tby, tbz,
714  rho_u_rhs, rho_v_rhs, rho_w_rhs,
715  tau11, tau22, tau33,
717  detJ_arr, stretched_dz_d, dxInv,
718  mf_mx, mf_ux, mf_vx,
719  mf_my, mf_uy, mf_vy,
720  l_use_stretched_dz,
721  l_use_terrain_fitted_coords);
722  } else {
723  DiffusionSrcForMom_EB(mfi, domain, tbx, tby, tbz,
724  rho_u_rhs, rho_v_rhs, rho_w_rhs,
725  u, v, w,
726  tau11, tau22, tau33,
727  tau12, tau13, tau23,
729  dx, dxInv,
730  mf_mx, mf_ux, mf_vx,
731  mf_my, mf_uy, mf_vy,
732  solverChoice, ebfact, bc_ptr_d);
733  }
734  }
735 
736  auto abl_pressure_grad = solverChoice.abl_pressure_grad;
737 
738  ParallelFor(tbx, tby,
739  [=] AMREX_GPU_DEVICE (int i, int j, int k)
740  { // x-momentum equation
741 
742  // Note that gradp arrays now carry the map factor in them
743 
744  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) : zero;
745 
746  rho_u_rhs(i, j, k) += (-gpx_arr(i,j,k) - abl_pressure_grad[0]) / (one + q) + xmom_src_arr(i,j,k);
747 
748  if (l_moving_terrain) {
749  Real h_zeta = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd);
750  rho_u_rhs(i, j, k) *= h_zeta;
751  }
752 
753  if ( l_anelastic && (nrk == 1) ) {
754  rho_u_rhs(i,j,k) *= myhalf;
755  rho_u_rhs(i,j,k) += myhalf / dt * (rho_u(i,j,k) - rho_u_old(i,j,k));
756  }
757  },
758  [=] AMREX_GPU_DEVICE (int i, int j, int k)
759  { // y-momentum equation
760 
761  // Note that gradp arrays now carry the map factor in them
762 
763  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) : zero;
764 
765  rho_v_rhs(i, j, k) += (-gpy_arr(i,j,k) - abl_pressure_grad[1]) / (one + q) + ymom_src_arr(i,j,k);
766 
767  if (l_moving_terrain) {
768  Real h_zeta = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd);
769  rho_v_rhs(i, j, k) *= h_zeta;
770  }
771 
772  if ( l_anelastic && (nrk == 1) ) {
773  rho_v_rhs(i,j,k) *= myhalf;
774  rho_v_rhs(i,j,k) += myhalf / dt * (rho_v(i,j,k) - rho_v_old(i,j,k));
775  }
776  });
777 
778  // *****************************************************************************
779  // Zero out source terms for x- and y- momenta if at walls or inflow
780  // We need to do this -- even though we call the boundary conditions later --
781  // because the slow source is used to update the state in the fast interpolater.
782  // *****************************************************************************
783  if (bx.smallEnd(0) == domain.smallEnd(0)) {
784  Box lo_x_dom_face(bx); lo_x_dom_face.setBig(0,bx.smallEnd(0));
785  if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir) {
786  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
787  rho_u_rhs(i,j,k) = zero;
788  });
789  } else if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir_upwind) {
790  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
791  if (u(i,j,k) >= zero) {
792  rho_u_rhs(i,j,k) = zero;
793  }
794  });
795  }
796  }
797  if (bx.bigEnd(0) == domain.bigEnd(0)) {
798  Box hi_x_dom_face(bx); hi_x_dom_face.setSmall(0,bx.bigEnd(0)+1); hi_x_dom_face.setBig(0,bx.bigEnd(0)+1);
799  if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir) {
800  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
801  rho_u_rhs(i,j,k) = zero;
802  });
803  } else if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir_upwind) {
804  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
805  if (u(i,j,k) <= zero) {
806  rho_u_rhs(i,j,k) = zero;
807  }
808  });
809  }
810  }
811  if (bx.smallEnd(1) == domain.smallEnd(1)) {
812  Box lo_y_dom_face(bx); lo_y_dom_face.setBig(1,bx.smallEnd(1));
813  if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir) {
814  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
815  rho_v_rhs(i,j,k) = zero;
816  });
817  } else if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir_upwind) {
818  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
819  if (v(i,j,k) >= zero) {
820  rho_v_rhs(i,j,k) = zero;
821  }
822  });
823  }
824  }
825  if (bx.bigEnd(1) == domain.bigEnd(1)) {
826  Box hi_y_dom_face(bx); hi_y_dom_face.setSmall(1,bx.bigEnd(1)+1); hi_y_dom_face.setBig(1,bx.bigEnd(1)+1);
827  if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir) {
828  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
829  rho_v_rhs(i,j,k) = zero;
830  });
831  } else if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir_upwind) {
832  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
833  if (v(i,j,k) <= zero) {
834  rho_v_rhs(i,j,k) = zero;
835  }
836  });
837  }
838  }
839 
840  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
841  { // z-momentum equation
842 
843  Real gpz = gpz_arr(i,j,k);
844 
845  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i,j,k-1)) : zero;
846 
847  rho_w_rhs(i, j, k) += (-gpz - abl_pressure_grad[2] + buoyancy_arr(i,j,k)) / (one + q) + zmom_src_arr(i,j,k);
848 
849  if (l_moving_terrain) {
850  rho_w_rhs(i, j, k) *= myhalf * (detJ_arr(i,j,k) + detJ_arr(i,j,k-1));
851  }
852  });
853 
854  auto const lo = lbound(bx);
855  auto const hi = ubound(bx);
856 
857  // Note: the logic below assumes no tiling in z!
858  if (level > 0) {
859 
860  const Array4<const Real>& rho_w_rhs_crse = zmom_crse_rhs->const_array(mfi);
861 
862  Box b2d = bx; b2d.setRange(2,0);
863 
864  if (lo.z > klo) {
865  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // bottom of box but not of domain
866  {
867  rho_w_rhs(i,j,lo.z) = rho_w_rhs_crse(i,j,lo.z);
868  });
869  }
870 
871  if (hi.z < khi+1) {
872  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // top of box but not of domain
873  {
874  rho_w_rhs(i,j,hi.z+1) = rho_w_rhs_crse(i,j,hi.z+1);
875  });
876  }
877  }
878 
879  {
880  BL_PROFILE("slow_rhs_pre_fluxreg");
881  // We only add to the flux registers in the final RK step
882  // NOTE: for now we are only refluxing density not (rho theta) since the latter seems to introduce
883  // a problem at top and bottom boundaries
884  if (l_reflux) {
885  int strt_comp_reflux = (l_fixed_rho) ? 1 : 0;
886  int num_comp_reflux = 1;
887  if (level < finest_level) {
888  fr_as_crse->CrseAdd(mfi,
889  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
890  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
891  }
892  if (level > 0) {
893  fr_as_fine->FineAdd(mfi,
894  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
895  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
896  }
897 
898  } // two-way coupling
899  } // end profile
900  } // mfi
901  } // OMP
902 }
void AdvectionSrcForRho(const amrex::Box &bx, const amrex::Array4< amrex::Real > &src, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &omega, const amrex::Array4< amrex::Real > &avg_xmom, const amrex::Array4< amrex::Real > &avg_ymom, const amrex::Array4< amrex::Real > &avg_zmom, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_arr, const bool fixed_rho)
void AdvectionSrcForScalars(const amrex::Box &bx, const int icomp, const int ncomp, const amrex::Array4< const amrex::Real > &avg_xmom, const amrex::Array4< const amrex::Real > &avg_ymom, const amrex::Array4< const amrex::Real > &avg_zmom, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &src, const amrex::Array4< const amrex::Real > &vf_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_my, const AdvType horiz_adv_type, const AdvType vert_adv_type, const amrex::Real horiz_upw_frac, const amrex::Real vert_upw_frac, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_arr, const amrex::Box &domain, const amrex::BCRec *bc_ptr_h)
void AdvectionSrcForMom(const amrex::MFIter &mfi, const amrex::Box &bx, const amrex::Box &bxx, const amrex::Box &bxy, const amrex::Box &bxz, const amrex::Vector< amrex::Box > &bxx_grown, const amrex::Vector< amrex::Box > &bxy_grown, const amrex::Vector< amrex::Box > &bxz_grown, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &rho, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &w, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &Omega, const amrex::Array4< const amrex::Real > &z_nd, const amrex::Array4< const amrex::Real > &ax, const amrex::Array4< const amrex::Real > &ay, const amrex::Array4< const amrex::Real > &az, const amrex::Array4< const amrex::Real > &detJ, amrex::Gpu::DeviceVector< amrex::Real > &stretched_dz_d, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, const AdvType horiz_adv_type, const AdvType vert_adv_type, const amrex::Real horiz_upw_frac, const amrex::Real vert_upw_frac, MeshType &mesh_type, TerrainType &terrain_type, const eb_ &ebfact, amrex::GpuArray< amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_u_arr, amrex::GpuArray< amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_v_arr, amrex::GpuArray< amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_w_arr, const amrex::Vector< amrex::iMultiFab > &physbnd_mask, const bool already_on_centroids, const int lo_z_face, const int hi_z_face, const amrex::Box &domain, const amrex::BCRec *bc_ptr_h)
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau32
Definition: ERF_DataStruct.H:32
@ tau31
Definition: ERF_DataStruct.H:32
@ tau21
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
@ nvars
Definition: ERF_DataStruct.H:98
@ v_x
Definition: ERF_DataStruct.H:24
@ u_y
Definition: ERF_DataStruct.H:25
@ v_y
Definition: ERF_DataStruct.H:25
@ m_y
Definition: ERF_DataStruct.H:25
@ u_x
Definition: ERF_DataStruct.H:24
@ m_x
Definition: ERF_DataStruct.H:24
void DiffusionSrcForMom(const amrex::Box &bxx, const amrex::Box &bxy, const amrex::Box &bxz, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &tau11, const amrex::Array4< const amrex::Real > &tau22, const amrex::Array4< const amrex::Real > &tau33, const amrex::Array4< const amrex::Real > &tau12, const amrex::Array4< const amrex::Real > &tau21, const amrex::Array4< const amrex::Real > &tau13, const amrex::Array4< const amrex::Real > &tau31, const amrex::Array4< const amrex::Real > &tau23, const amrex::Array4< const amrex::Real > &tau32, const amrex::Array4< const amrex::Real > &detJ_arr, const amrex::Gpu::DeviceVector< amrex::Real > &stretched_dz_d, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, const bool use_stretched_dz, const bool use_variable_dz)
void DiffusionSrcForState_S(const amrex::Box &bx, const amrex::Box &domain, int start_comp, int num_comp, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< amrex::Real > &xflux, const amrex::Array4< amrex::Real > &yflux, const amrex::Array4< amrex::Real > &zflux, const amrex::Gpu::DeviceVector< amrex::Real > &stretched_dz_d, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &SmnSmn_a, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, amrex::Array4< amrex::Real > &hfx_z, amrex::Array4< amrex::Real > &qfx1_z, amrex::Array4< amrex::Real > &qfx2_z, amrex::Array4< amrex::Real > &diss, const amrex::Array4< const amrex::Real > &mu_turb, const SolverChoice &solverChoice, const int level, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > grav_gpu, const amrex::BCRec *bc_ptr, const bool use_SurfLayer, const amrex::Real implicit_fac)
void DiffusionSrcForState_T(const amrex::Box &bx, const amrex::Box &domain, int start_comp, int num_comp, const bool &rotate, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< amrex::Real > &xflux, const amrex::Array4< amrex::Real > &yflux, const amrex::Array4< amrex::Real > &zflux, const amrex::Array4< const amrex::Real > &z_nd, const amrex::Array4< const amrex::Real > &z_cc, const amrex::Array4< const amrex::Real > &ax, const amrex::Array4< const amrex::Real > &ay, const amrex::Array4< const amrex::Real > &az, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &SmnSmn_a, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, amrex::Array4< amrex::Real > &hfx_x, amrex::Array4< amrex::Real > &hfx_y, amrex::Array4< amrex::Real > &hfx_z, amrex::Array4< amrex::Real > &qfx1_x, amrex::Array4< amrex::Real > &qfx1_y, amrex::Array4< amrex::Real > &qfx1_z, amrex::Array4< amrex::Real > &qfx2_z, amrex::Array4< amrex::Real > &diss, const amrex::Array4< const amrex::Real > &mu_turb, const SolverChoice &solverChoice, const int level, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > grav_gpu, const amrex::BCRec *bc_ptr, const bool use_SurfLayer, const amrex::Real implicit_fac)
void DiffusionSrcForState_N(const amrex::Box &bx, const amrex::Box &domain, int start_comp, int num_comp, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< amrex::Real > &xflux, const amrex::Array4< amrex::Real > &yflux, const amrex::Array4< amrex::Real > &zflux, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &SmnSmn_a, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, amrex::Array4< amrex::Real > &hfx_z, amrex::Array4< amrex::Real > &qfx1_z, amrex::Array4< amrex::Real > &qfx2_z, amrex::Array4< amrex::Real > &diss, const amrex::Array4< const amrex::Real > &mu_turb, const SolverChoice &solverChoice, const int level, const amrex::Array4< const amrex::Real > &tm_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > grav_gpu, const amrex::BCRec *bc_ptr, const bool use_SurfLayer, const amrex::Real implicit_fac)
void DiffusionSrcForMom_EB(const amrex::MFIter &mfi, [[maybe_unused]] const amrex::Box &domain, const amrex::Box &bxx, const amrex::Box &bxy, const amrex::Box &bxz, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &w_arr, const amrex::Array4< const amrex::Real > &tau11, const amrex::Array4< const amrex::Real > &tau22, const amrex::Array4< const amrex::Real > &tau33, const amrex::Array4< const amrex::Real > &tau12, const amrex::Array4< const amrex::Real > &tau13, const amrex::Array4< const amrex::Real > &tau23, const amrex::Array4< const amrex::Real > &tau_eb13, const amrex::Array4< const amrex::Real > &tau_eb23, const amrex::Real *dx_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_ux, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vy, const SolverChoice &solverChoice, const eb_ &ebfact, [[maybe_unused]] const amrex::BCRec *bc_ptr)
void DiffusionSrcForState_EB(const amrex::Box &bx, const amrex::Box &domain, int start_comp, int num_comp, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< amrex::Real > &xflux, const amrex::Array4< amrex::Real > &yflux, const amrex::Array4< amrex::Real > &zflux, const amrex::Array4< const amrex::EBCellFlag > &cfg_arr, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Array4< const amrex::Real > &detJ, const amrex::Array4< const amrex::Real > &barea_arr, const amrex::Array4< const amrex::Real > &bcent_arr, const amrex::Real *dx_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, amrex::Array4< amrex::Real > &hfx_z, amrex::Array4< amrex::Real > &qfx1_z, amrex::Array4< amrex::Real > &qfx2_z, amrex::Array4< amrex::Real > &hfx_EB, const amrex::Array4< const amrex::Real > &mu_turb, const SolverChoice &solverChoice, const int level, const amrex::BCRec *bc_ptr, const bool use_SurfLayer)
void EBAdvectionSrcForScalars(const amrex::Box &bx, const int icomp, const int ncomp, const amrex::Array4< const amrex::Real > &avg_xmom, const amrex::Array4< const amrex::Real > &avg_ymom, const amrex::Array4< const amrex::Real > &avg_zmom, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &src, const amrex::Array4< const int > &mask_arr, const amrex::Array4< const amrex::EBCellFlag > &cfg_arr, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Array4< const amrex::Real > &fcx_arr, const amrex::Array4< const amrex::Real > &fcy_arr, const amrex::Array4< const amrex::Real > &fcz_arr, const amrex::Array4< const amrex::Real > &vf_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_my, const AdvType horiz_adv_type, const AdvType vert_adv_type, const amrex::Real horiz_upw_frac, const amrex::Real vert_upw_frac, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_arr, const amrex::Box &domain, const amrex::BCRec *bc_ptr_h, bool already_on_centroids)
void EBAdvectionSrcForRho(const amrex::Box &bx, const amrex::Array4< amrex::Real > &src, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &omega, const amrex::Array4< amrex::Real > &avg_xmom, const amrex::Array4< amrex::Real > &avg_ymom, const amrex::Array4< amrex::Real > &avg_zmom, const amrex::Array4< const int > &mask_arr, const amrex::Array4< const amrex::EBCellFlag > &cfg_arr, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Array4< const amrex::Real > &fcx_arr, const amrex::Array4< const amrex::Real > &fcy_arr, const amrex::Array4< const amrex::Real > &fcz_arr, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &mf_mx, const amrex::Array4< const amrex::Real > &mf_my, const amrex::Array4< const amrex::Real > &mf_uy, const amrex::Array4< const amrex::Real > &mf_vx, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_arr, const bool fixed_rho, bool already_on_centroids)
@ tau_eb23
Definition: ERF_EBStruct.H:16
@ tau_eb13
Definition: ERF_EBStruct.H:16
const Real l_vert_implicit_fac
Definition: ERF_ImplicitPost.H:5
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
AdvType
Definition: ERF_IndexDefines.H:239
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
void erf_make_tau_terms(int level, int nrk, const Vector< BCRec > &domain_bcs_type_h, const MultiFab &z_phys_nd, Vector< MultiFab > &S_data, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &zvel, Vector< std::unique_ptr< MultiFab >> &Tau_lev, Vector< std::unique_ptr< MultiFab >> &Tau_corr_lev, MultiFab *SmnSmn, MultiFab *eddyDiffs, const Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &, Gpu::DeviceVector< Real > &stretched_dz_d, const MultiFab &detJ, Vector< std::unique_ptr< MultiFab >> &mapfac, const MultiFab &ax, const MultiFab &ay, const MultiFab &az, const eb_ &ebfact)
Definition: ERF_MakeTauTerms.cpp:12
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real OmegaFromW(int &i, int &j, int &k, amrex::Real w, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:414
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:104
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:144
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
const std::unique_ptr< amrex::EBFArrayBoxFactory > & get_const_factory() const noexcept
Definition: ERF_EB.H:46
@ yvel_bc
Definition: ERF_IndexDefines.H:88
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ ext_dir
Definition: ERF_IndexDefines.H:227
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:235
@ gpz
Definition: ERF_IndexDefines.H:170
@ gpy
Definition: ERF_IndexDefines.H:169
@ gpx
Definition: ERF_IndexDefines.H:168
@ NumTypes
Definition: ERF_IndexDefines.H:180
@ ymom
Definition: ERF_IndexDefines.H:178
@ cons
Definition: ERF_IndexDefines.H:176
@ zmom
Definition: ERF_IndexDefines.H:179
@ xmom
Definition: ERF_IndexDefines.H:177
@ qt
Definition: ERF_Kessler.H:27
@ ng
Definition: ERF_Morrison.H:48
@ xvel
Definition: ERF_IndexDefines.H:159
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160
AdvType dycore_vert_adv_type
Definition: ERF_AdvStruct.H:420
amrex::Real dycore_vert_upw_frac
Definition: ERF_AdvStruct.H:430
AdvType dycore_horiz_adv_type
Definition: ERF_AdvStruct.H:419
amrex::Real dycore_horiz_upw_frac
Definition: ERF_AdvStruct.H:429
Definition: ERF_DiffStruct.H:19
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
static MeshType mesh_type
Definition: ERF_DataStruct.H:1079
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1088
amrex::Real gravity
Definition: ERF_DataStruct.H:1154
bool use_shoc
Definition: ERF_DataStruct.H:1187
amrex::Vector< amrex::Real > vert_implicit_fac
Definition: ERF_DataStruct.H:1113
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > abl_pressure_grad
Definition: ERF_DataStruct.H:1213
bool implicit_thermal_diffusion
Definition: ERF_DataStruct.H:1115
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1091
amrex::Vector< int > anelastic
Definition: ERF_DataStruct.H:1097
AdvChoice advChoice
Definition: ERF_DataStruct.H:1087
MoistureType moisture_type
Definition: ERF_DataStruct.H:1206
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1067
amrex::Vector< int > fixed_density
Definition: ERF_DataStruct.H:1098
bool use_rotate_surface_flux
Definition: ERF_DataStruct.H:1184
CouplingType coupling_type
Definition: ERF_DataStruct.H:1205
Definition: ERF_TurbStruct.H:42
PBLType pbl_type
Definition: ERF_TurbStruct.H:418
bool use_keqn
Definition: ERF_TurbStruct.H:425
RANSType rans_type
Definition: ERF_TurbStruct.H:413
LESType les_type
Definition: ERF_TurbStruct.H:371
bool use_kturb
Definition: ERF_TurbStruct.H:424
Here is the call graph for this function: