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< 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< 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] : zero;
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  }
247  }
248 #endif
249  } // l_use_diff
250 
251  // This is just cautionary to deal with grid boundaries that aren't domain boundaries
252  S_rhs[IntVars::zmom].setVal(0);
253 
254  // *****************************************************************************
255  // Define updates and fluxes in the current RK stage
256  // *****************************************************************************
257  // Cell-centered masks for EB (used for flux interpolation)
258  bool already_on_centroids = false;
259  Vector<iMultiFab> physbnd_mask;
260  physbnd_mask.resize(IntVars::NumTypes);
261  if (l_use_eb) {
262  physbnd_mask[IntVars::cons].define(S_data[IntVars::cons].boxArray(), S_data[IntVars::cons].DistributionMap(), 1, 1);
263  physbnd_mask[IntVars::cons].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
264  // physbnd_mask[IntVars::cons].FillBoundary(geom.periodicity());
265  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
266  physbnd_mask[1+dir].define(S_data[1+dir].boxArray(), S_data[1+dir].DistributionMap(), 1, 1);
267  physbnd_mask[1+dir].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
268  // physbnd_mask[1+dir].FillBoundary(geom.periodicity());
269  }
270  }
271 
272 #ifdef _OPENMP
273 #pragma omp parallel if (Gpu::notInLaunchRegion())
274 #endif
275  {
276  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
277  {
278  Box bx = mfi.tilebox();
279  Box tbx = mfi.nodaltilebox(0);
280  Box tby = mfi.nodaltilebox(1);
281  Box tbz = mfi.nodaltilebox(2);
282 
283  // Boxes for momentum fluxes
284  Vector<Box> tbx_grown(AMREX_SPACEDIM);
285  Vector<Box> tby_grown(AMREX_SPACEDIM);
286  Vector<Box> tbz_grown(AMREX_SPACEDIM);
287  if (l_use_eb) {
288  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
289  tbx_grown[dir] = tbx;
290  tby_grown[dir] = tby;
291  tbz_grown[dir] = tbz;
292  IntVect iv(1, 1, 1);
293  iv[dir] = 0;
294  tbx_grown[dir] = (tbx_grown[dir].growHi(dir,1)).grow(iv);
295  tby_grown[dir] = (tby_grown[dir].growHi(dir,1)).grow(iv);
296  tbz_grown[dir] = (tbz_grown[dir].growHi(dir,1)).grow(iv);
297  }
298  }
299 
300  // We don't compute a source term for z-momentum on the bottom or top domain boundary
301  if (tbz.smallEnd(2) == domain.smallEnd(2)) {
302  tbz.growLo(2,-1);
303  }
304  if (tbz.bigEnd(2) == domain.bigEnd(2)+1) {
305  tbz.growHi(2,-1);
306  }
307 
308  const Array4<const Real> & cell_data = S_data[IntVars::cons].array(mfi);
309  const Array4<const Real> & cell_prim = S_prim.array(mfi);
310  const Array4<Real> & cell_rhs = S_rhs[IntVars::cons].array(mfi);
311 
312  const Array4<const Real> & cell_old = S_old[IntVars::cons].array(mfi);
313 
314  const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
315  const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
316  const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
317  const Array4<Real const>& buoyancy_arr = buoyancy.const_array(mfi);
318 
319  const Array4<Real const>& gpx_arr = gradp[GpVars::gpx].const_array(mfi);
320  const Array4<Real const>& gpy_arr = gradp[GpVars::gpy].const_array(mfi);
321  const Array4<Real const>& gpz_arr = gradp[GpVars::gpz].const_array(mfi);
322 
323  const Array4<Real const>& qt_arr = qt.const_array(mfi);
324 
325  const Array4<Real>& rho_u_old = S_old[IntVars::xmom].array(mfi);
326  const Array4<Real>& rho_v_old = S_old[IntVars::ymom].array(mfi);
327 
328  if (l_anelastic) {
329  // When anelastic we must reset these to 0 each RK step
330  avg_xmom[mfi].template setVal<RunOn::Device>(0,tbx);
331  avg_ymom[mfi].template setVal<RunOn::Device>(0,tby);
332  avg_zmom[mfi].template setVal<RunOn::Device>(0,tbz);
333  }
334 
335  Array4<Real> avg_xmom_arr = avg_xmom.array(mfi);
336  Array4<Real> avg_ymom_arr = avg_ymom.array(mfi);
337  Array4<Real> avg_zmom_arr = avg_zmom.array(mfi);
338 
339  const Array4<const Real> & u = xvel.array(mfi);
340  const Array4<const Real> & v = yvel.array(mfi);
341  const Array4<const Real> & w = zvel.array(mfi);
342 
343  const Array4<const Real>& rho_u = S_data[IntVars::xmom].array(mfi);
344  const Array4<const Real>& rho_v = S_data[IntVars::ymom].array(mfi);
345  const Array4<const Real>& rho_w = S_data[IntVars::zmom].array(mfi);
346 
347  // Map factors
348  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
349  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
350  const Array4<const Real>& mf_vx = mapfac[MapFacType::v_x]->const_array(mfi);
351  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
352  const Array4<const Real>& mf_uy = mapfac[MapFacType::u_y]->const_array(mfi);
353  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
354 
355  const Array4< Real>& omega_arr = Omega.array(mfi);
356 
357  Array4<const Real> z_t;
358  if (z_t_mf) {
359  z_t = z_t_mf->array(mfi);
360  } else {
361  z_t = Array4<const Real>{};
362  }
363 
364  const Array4<Real>& rho_u_rhs = S_rhs[IntVars::xmom].array(mfi);
365  const Array4<Real>& rho_v_rhs = S_rhs[IntVars::ymom].array(mfi);
366  const Array4<Real>& rho_w_rhs = S_rhs[IntVars::zmom].array(mfi);
367 
368  const Array4<Real const>& mu_turb = l_use_turb ? eddyDiffs->const_array(mfi) : Array4<const Real>{};
369 
370  // Terrain metrics
371  const Array4<const Real>& z_nd = z_phys_nd.const_array(mfi);
372  const Array4<const Real>& z_cc = z_phys_cc.const_array(mfi);
373 
374  // *****************************************************************************
375  // Define flux arrays for use in advection
376  // *****************************************************************************
377  std::array<FArrayBox,AMREX_SPACEDIM> flux;
378  std::array<FArrayBox,AMREX_SPACEDIM> flux_u;
379  std::array<FArrayBox,AMREX_SPACEDIM> flux_v;
380  std::array<FArrayBox,AMREX_SPACEDIM> flux_w;
381 
382  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
383  if (!l_use_eb) {
384  flux[dir].resize(surroundingNodes(bx,dir),2,The_Async_Arena());
385  } else {
386  flux[dir].resize(surroundingNodes(bx,dir).grow(1),2,The_Async_Arena());
387  }
388  flux[dir].setVal<RunOn::Device>(0);
389  }
390  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
391  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
392 
393  // Define flux arrays for momentum variables (used only for EB now)
394  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_u_arr{};
395  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_v_arr{};
396  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_w_arr{};
397 
398  if (l_use_eb) {
399  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
400  flux_u[dir].resize(tbx_grown[dir],1,The_Async_Arena());
401  flux_v[dir].resize(tby_grown[dir],1,The_Async_Arena());
402  flux_w[dir].resize(tbz_grown[dir],1,The_Async_Arena());
403  flux_u[dir].setVal<RunOn::Device>(0);
404  flux_v[dir].setVal<RunOn::Device>(0);
405  flux_w[dir].setVal<RunOn::Device>(0);
406  flx_u_arr[dir] = flux_u[dir].array();
407  flx_v_arr[dir] = flux_v[dir].array();
408  flx_w_arr[dir] = flux_w[dir].array();
409  }
410  }
411 
412  // *****************************************************************************
413  // Contravariant flux field
414  // *****************************************************************************
415  {
416  BL_PROFILE("slow_rhs_making_omega");
417  IntVect nGrowVect = (l_use_eb)
418  ? IntVect(AMREX_D_DECL(2, 2, 2)) : IntVect(AMREX_D_DECL(1, 1, 1));
419  Box gbxo = surroundingNodes(bx,2); gbxo.grow(nGrowVect);
420  //
421  // Now create Omega with momentum (not velocity) with z_t subtracted if moving terrain
422  // ONLY if not doing anelastic + terrain -- in that case Omega will be defined coming
423  // out of the projection
424  //
425  if (!l_use_terrain_fitted_coords) {
426  ParallelFor(gbxo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
427  omega_arr(i,j,k) = rho_w(i,j,k);
428  });
429 
430  } else {
431 
432  Box gbxo_lo = gbxo; gbxo_lo.setBig(2,domain.smallEnd(2));
433  int lo_z_face = domain.smallEnd(2);
434  if (gbxo_lo.smallEnd(2) <= lo_z_face) {
435  ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
436  omega_arr(i,j,k) = zero;
437  });
438  }
439  Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
440  int hi_z_face = domain.bigEnd(2)+1;
441  if (gbxo_hi.bigEnd(2) >= hi_z_face) {
442  ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
443  omega_arr(i,j,k) = rho_w(i,j,k);
444  });
445  }
446 
447  if (z_t) { // Note we never do anelastic with moving terrain
448  Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
449  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
450  // We define rho on the z-face the same way as in MomentumToVelocity/VelocityToMomentum
451  Real rho_at_face = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp));
452  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
453  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv) -
454  rho_at_face * z_t(i,j,k);
455  });
456  } else {
457  Box gbxo_mid = gbxo;
458  if (gbxo_mid.smallEnd(2) <= domain.smallEnd(2)) {
459  gbxo_mid.setSmall(2,1);
460  }
461  if (gbxo_mid.bigEnd(2) >= domain.bigEnd(2)+1) {
462  gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
463  }
464  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
465  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
466  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv);
467  });
468  }
469  }
470  } // end profile
471 
472  // *****************************************************************************
473  // Diffusive terms (pre-computed above)
474  // *****************************************************************************
475  // No terrain diffusion
476  Array4<Real> tau11,tau22,tau33;
477  Array4<Real> tau12,tau13,tau23;
478  if (Tau_lev[TauType::tau11]) {
479  tau11 = Tau_lev[TauType::tau11]->array(mfi); tau22 = Tau_lev[TauType::tau22]->array(mfi);
480  tau33 = Tau_lev[TauType::tau33]->array(mfi); tau12 = Tau_lev[TauType::tau12]->array(mfi);
481  tau13 = Tau_lev[TauType::tau13]->array(mfi); tau23 = Tau_lev[TauType::tau23]->array(mfi);
482  } else {
483  tau11 = Array4<Real>{}; tau22 = Array4<Real>{}; tau33 = Array4<Real>{};
484  tau12 = Array4<Real>{}; tau13 = Array4<Real>{}; tau23 = Array4<Real>{};
485  }
486  // Terrain diffusion
487  Array4<Real> tau21,tau31,tau32;
488  if (Tau_lev[TauType::tau21]) {
489  tau21 = Tau_lev[TauType::tau21]->array(mfi);
490  tau31 = Tau_lev[TauType::tau31]->array(mfi);
491  tau32 = Tau_lev[TauType::tau32]->array(mfi);
492  } else {
493  tau21 = Array4<Real>{}; tau31 = Array4<Real>{}; tau32 = Array4<Real>{};
494  }
495 
496  // EB surface layer fluxes
497  Array4<Real> u_tau_eb13, u_tau_eb23;
498  Array4<Real> v_tau_eb13, v_tau_eb23;
499  Array4<Real> w_tau_eb13, w_tau_eb23;
500  if (l_use_eb) {
501  EBChoice ebChoice = solverChoice.ebChoice;
503  u_tau_eb13 = Tau_EB[EBTauType::tau_eb13][EBGridType::xface]->array(mfi);
504  u_tau_eb23 = Tau_EB[EBTauType::tau_eb23][EBGridType::xface]->array(mfi);
505  v_tau_eb13 = Tau_EB[EBTauType::tau_eb13][EBGridType::yface]->array(mfi);
506  v_tau_eb23 = Tau_EB[EBTauType::tau_eb23][EBGridType::yface]->array(mfi);
507  w_tau_eb13 = Tau_EB[EBTauType::tau_eb13][EBGridType::zface]->array(mfi);
508  w_tau_eb23 = Tau_EB[EBTauType::tau_eb23][EBGridType::zface]->array(mfi);
509  }
510  }
511 
512  // Strain magnitude
513  Array4<Real> SmnSmn_a;
514  if (l_need_SmnSmn) {
515  SmnSmn_a = SmnSmn->array(mfi);
516  } else {
517  SmnSmn_a = Array4<Real>{};
518  }
519 
520  // *****************************************************************************
521  // Define updates in the RHS of continuity and potential temperature equations
522  // *****************************************************************************
523  bool l_eb_terrain_cc = false; // EB terrain on cell-centered grid
524  Array4<const int> mask_arr{};
525  Array4<const EBCellFlag> cfg_arr{};
526  Array4<const Real> ax_arr{};
527  Array4<const Real> ay_arr{};
528  Array4<const Real> az_arr{};
529  Array4<const Real> fcx_arr{};
530  Array4<const Real> fcy_arr{};
531  Array4<const Real> fcz_arr{};
532  Array4<const Real> detJ_arr{};
533  Array4<const Real> barea_arr{};
534  Array4<const Real> bcent_arr{};
535 
536  if (l_use_eb)
537  {
538  EBCellFlagFab const& cfg = (ebfact.get_const_factory())->getMultiEBCellFlagFab()[mfi];
539  cfg_arr = cfg.const_array();
540  if (cfg.getType(bx) == FabType::singlevalued) {
541  l_eb_terrain_cc = true;
542  ax_arr = (ebfact.get_const_factory())->getAreaFrac()[0]->const_array(mfi);
543  ay_arr = (ebfact.get_const_factory())->getAreaFrac()[1]->const_array(mfi);
544  az_arr = (ebfact.get_const_factory())->getAreaFrac()[2]->const_array(mfi);
545  fcx_arr = (ebfact.get_const_factory())->getFaceCent()[0]->const_array(mfi);
546  fcy_arr = (ebfact.get_const_factory())->getFaceCent()[1]->const_array(mfi);
547  fcz_arr = (ebfact.get_const_factory())->getFaceCent()[2]->const_array(mfi);
548  detJ_arr = (ebfact.get_const_factory())->getVolFrac().const_array(mfi);
549  mask_arr = physbnd_mask[IntVars::cons].const_array(mfi);
550  barea_arr = (ebfact.get_const_factory())->getBndryArea().const_array(mfi);
551  bcent_arr = (ebfact.get_const_factory())->getBndryCent().const_array(mfi);
552  } else {
553  ax_arr = ax.const_array(mfi);
554  ay_arr = ay.const_array(mfi);
555  az_arr = az.const_array(mfi);
556  detJ_arr = detJ.const_array(mfi);
557  }
558  } else {
559  ax_arr = ax.const_array(mfi);
560  ay_arr = ay.const_array(mfi);
561  az_arr = az.const_array(mfi);
562  detJ_arr = detJ.const_array(mfi);
563  }
564 
565  int icomp = RhoTheta_comp; int ncomp = 1;
566  if (!l_eb_terrain_cc){
567  AdvectionSrcForRho( bx, cell_rhs,
568  rho_u, rho_v, omega_arr, // these are being used to build the fluxes
569  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr, // these are being defined from the fluxes
570  ax_arr, ay_arr, az_arr, detJ_arr,
571  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
572  flx_arr, l_fixed_rho);
573  AdvectionSrcForScalars(bx, icomp, ncomp,
574  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
575  cell_prim, cell_rhs,
576  detJ_arr, dxInv, mf_mx, mf_my,
577  l_horiz_adv_type, l_vert_adv_type,
578  l_horiz_upw_frac, l_vert_upw_frac,
579  flx_arr, domain, bc_ptr_h);
580  } else {
581  EBAdvectionSrcForRho(bx, cell_rhs,
582  rho_u, rho_v, omega_arr,
583  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
584  mask_arr, cfg_arr,
585  ax_arr, ay_arr, az_arr,
586  fcx_arr, fcy_arr, fcz_arr, detJ_arr,
587  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
588  flx_arr, l_fixed_rho,
589  already_on_centroids);
590  EBAdvectionSrcForScalars(bx, icomp, ncomp,
591  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
592  cell_prim, cell_rhs,
593  mask_arr, cfg_arr, ax_arr, ay_arr, az_arr,
594  fcx_arr, fcy_arr, fcz_arr,
595  detJ_arr, dxInv, mf_mx, mf_my,
596  l_horiz_adv_type, l_vert_adv_type,
597  l_horiz_upw_frac, l_vert_upw_frac,
598  flx_arr, domain, bc_ptr_h,
599  already_on_centroids);
600  }
601 
602  if (l_use_diff) {
603  Array4<Real> diffflux_x = dflux_x->array(mfi);
604  Array4<Real> diffflux_y = dflux_y->array(mfi);
605  Array4<Real> diffflux_z = dflux_z->array(mfi);
606 
607  Array4<Real> hfx_x = Hfx1->array(mfi);
608  Array4<Real> hfx_y = Hfx2->array(mfi);
609  Array4<Real> hfx_z = Hfx3->array(mfi);
610  Array4<Real> hfx_EB{};
611  if (l_use_eb) {
612  hfx_EB = Hfx3_EB->array(mfi);
613  }
614 
615  Array4<Real> q1fx_x = (Q1fx1) ? Q1fx1->array(mfi) : Array4<Real>{};
616  Array4<Real> q1fx_y = (Q1fx2) ? Q1fx2->array(mfi) : Array4<Real>{};
617  Array4<Real> q1fx_z = (Q1fx3) ? Q1fx3->array(mfi) : Array4<Real>{};
618 
619  Array4<Real> q2fx_z = (Q2fx3) ? Q2fx3->array(mfi) : Array4<Real>{};
620  Array4<Real> diss = Diss->array(mfi);
621 
622  const Array4<const Real> tm_arr = t_mean_mf ? t_mean_mf->const_array(mfi) : Array4<const Real>{};
623 
624  // NOTE: No diffusion for continuity, so n starts at one
625  int n_start = RhoTheta_comp;
626  int n_comp = 1;
627 
628  // For l_vert_implicit_fac > 0, we scale the rho*theta contribution
629  // by (1 - implicit_fac) and add in the implicit contribution with
630  // ERF_Implicit.H
631  if (l_use_stretched_dz) {
632  DiffusionSrcForState_S(bx, domain, n_start, n_comp, u, v,
633  cell_data, cell_prim, cell_rhs,
634  diffflux_x, diffflux_y, diffflux_z,
635  stretched_dz_d, dxInv, SmnSmn_a,
636  mf_mx, mf_ux, mf_vx,
637  mf_my, mf_uy, mf_vy,
638  hfx_z, q1fx_z, q2fx_z, diss,
639  mu_turb, solverChoice, level,
640  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
641  } else if (l_use_terrain_fitted_coords) {
642  DiffusionSrcForState_T(bx, domain, n_start, n_comp, l_rotate, u, v,
643  cell_data, cell_prim, cell_rhs,
644  diffflux_x, diffflux_y, diffflux_z,
645  z_nd, z_cc, ax_arr, ay_arr, az_arr, detJ_arr,
646  dxInv, SmnSmn_a,
647  mf_mx, mf_ux, mf_vx,
648  mf_my, mf_uy, mf_vy,
649  hfx_x, hfx_y, hfx_z, q1fx_x, q1fx_y, q1fx_z, q2fx_z, diss,
650  mu_turb, solverChoice, level,
651  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
652  } else if (l_use_eb) {
653  DiffusionSrcForState_EB(bx, domain, n_start, n_comp, u, v,
654  cell_data, cell_prim, cell_rhs,
655  diffflux_x, diffflux_y, diffflux_z,
656  cfg_arr, ax_arr, ay_arr, az_arr, detJ_arr,
657  barea_arr, bcent_arr,
658  dx, dxInv,
659  hfx_z, q1fx_z, q2fx_z, hfx_EB,
660  mu_turb, solverChoice, level,
661  bc_ptr_d, l_use_SurfLayer);
662  } else {
663  DiffusionSrcForState_N(bx, domain, n_start, n_comp, u, v,
664  cell_data, cell_prim, cell_rhs,
665  diffflux_x, diffflux_y, diffflux_z,
666  dxInv, SmnSmn_a,
667  mf_mx, mf_ux, mf_vx,
668  mf_my, mf_uy, mf_vy,
669  hfx_z, q1fx_z, q2fx_z, diss,
670  mu_turb, solverChoice, level,
671  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer, l_vert_implicit_fac);
672  }
673  }
674 
675  const Array4<Real const>& source_arr = cc_src.const_array(mfi);
676  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
677  {
678  cell_rhs(i,j,k,Rho_comp) += source_arr(i,j,k,Rho_comp);
679  cell_rhs(i,j,k,RhoTheta_comp) += source_arr(i,j,k,RhoTheta_comp);
680  });
681 
682  // If anelastic and in second RK stage, take average of old-time and new-time source
683  if ( l_anelastic && (nrk == 1) )
684  {
685  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
686  {
687  cell_rhs(i,j,k, Rho_comp) *= myhalf;
688  cell_rhs(i,j,k,RhoTheta_comp) *= myhalf;
689 
690  cell_rhs(i,j,k, Rho_comp) += myhalf / dt * (cell_data(i,j,k, Rho_comp) - cell_old(i,j,k, Rho_comp));
691  cell_rhs(i,j,k,RhoTheta_comp) += myhalf / dt * (cell_data(i,j,k,RhoTheta_comp) - cell_old(i,j,k,RhoTheta_comp));
692  });
693  }
694 
695  // *****************************************************************************
696  // Define updates in the RHS of {x, y, z}-momentum equations
697  // *****************************************************************************
698  int lo_z_face = domain.smallEnd(2);
699  int hi_z_face = domain.bigEnd(2)+1;
700 
701  AdvectionSrcForMom(mfi, bx, tbx, tby, tbz, tbx_grown, tby_grown, tbz_grown,
702  rho_u_rhs, rho_v_rhs, rho_w_rhs,
703  cell_data, u, v, w,
704  rho_u, rho_v, omega_arr,
705  z_nd, ax_arr, ay_arr, az_arr,
706  detJ_arr, stretched_dz_d,
707  dxInv, mf_mx, mf_ux, mf_vx, mf_my, mf_uy, mf_vy,
708  l_horiz_adv_type, l_vert_adv_type,
709  l_horiz_upw_frac, l_vert_upw_frac,
710  solverChoice.mesh_type, solverChoice.terrain_type,
711  ebfact, flx_u_arr, flx_v_arr, flx_w_arr,
712  physbnd_mask, already_on_centroids,
713  lo_z_face, hi_z_face, domain, bc_ptr_h);
714 
715  if (l_use_diff) {
716  // Note: tau** were calculated with calls to
717  // ComputeStress[Cons|Var]Visc_[N|S|T] in which ConsVisc ("constant
718  // viscosity") means that there is no contribution from a
719  // turbulence model. However, whether this field truly is constant
720  // depends on whether MolecDiffType is Constant or ConstantAlpha.
721  if (!l_use_eb) {
722  DiffusionSrcForMom(tbx, tby, tbz,
723  rho_u_rhs, rho_v_rhs, rho_w_rhs,
724  tau11, tau22, tau33,
726  detJ_arr, stretched_dz_d, dxInv,
727  mf_mx, mf_ux, mf_vx,
728  mf_my, mf_uy, mf_vy,
729  l_use_stretched_dz,
730  l_use_terrain_fitted_coords);
731  } else {
732  DiffusionSrcForMom_EB(mfi, domain, tbx, tby, tbz,
733  rho_u_rhs, rho_v_rhs, rho_w_rhs,
734  u, v, w,
735  tau11, tau22, tau33,
736  tau12, tau13, tau23,
737  u_tau_eb13, u_tau_eb23, v_tau_eb13, v_tau_eb23, w_tau_eb13, w_tau_eb23,
738  dx, dxInv,
739  mf_mx, mf_ux, mf_vx,
740  mf_my, mf_uy, mf_vy,
741  solverChoice, ebfact, bc_ptr_d);
742  }
743  }
744 
745  auto abl_pressure_grad = solverChoice.abl_pressure_grad;
746 
747  ParallelFor(tbx, tby,
748  [=] AMREX_GPU_DEVICE (int i, int j, int k)
749  { // x-momentum equation
750 
751  // Note that gradp arrays now carry the map factor in them
752 
753  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) : zero;
754 
755  rho_u_rhs(i, j, k) += (-gpx_arr(i,j,k) - abl_pressure_grad[0]) / (one + q) + xmom_src_arr(i,j,k);
756 
757  if (l_moving_terrain) {
758  Real h_zeta = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd);
759  rho_u_rhs(i, j, k) *= h_zeta;
760  }
761 
762  if ( l_anelastic && (nrk == 1) ) {
763  rho_u_rhs(i,j,k) *= myhalf;
764  rho_u_rhs(i,j,k) += myhalf / dt * (rho_u(i,j,k) - rho_u_old(i,j,k));
765  }
766  },
767  [=] AMREX_GPU_DEVICE (int i, int j, int k)
768  { // y-momentum equation
769 
770  // Note that gradp arrays now carry the map factor in them
771 
772  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) : zero;
773 
774  rho_v_rhs(i, j, k) += (-gpy_arr(i,j,k) - abl_pressure_grad[1]) / (one + q) + ymom_src_arr(i,j,k);
775 
776  if (l_moving_terrain) {
777  Real h_zeta = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd);
778  rho_v_rhs(i, j, k) *= h_zeta;
779  }
780 
781  if ( l_anelastic && (nrk == 1) ) {
782  rho_v_rhs(i,j,k) *= myhalf;
783  rho_v_rhs(i,j,k) += myhalf / dt * (rho_v(i,j,k) - rho_v_old(i,j,k));
784  }
785  });
786 
787  // *****************************************************************************
788  // Zero out source terms for x- and y- momenta if at walls or inflow
789  // We need to do this -- even though we call the boundary conditions later --
790  // because the slow source is used to update the state in the fast interpolater.
791  // *****************************************************************************
792  if (bx.smallEnd(0) == domain.smallEnd(0)) {
793  Box lo_x_dom_face(bx); lo_x_dom_face.setBig(0,bx.smallEnd(0));
794  if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir) {
795  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
796  rho_u_rhs(i,j,k) = zero;
797  });
798  } else if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir_upwind) {
799  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
800  if (u(i,j,k) >= zero) {
801  rho_u_rhs(i,j,k) = zero;
802  }
803  });
804  }
805  }
806  if (bx.bigEnd(0) == domain.bigEnd(0)) {
807  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);
808  if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir) {
809  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
810  rho_u_rhs(i,j,k) = zero;
811  });
812  } else if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir_upwind) {
813  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
814  if (u(i,j,k) <= zero) {
815  rho_u_rhs(i,j,k) = zero;
816  }
817  });
818  }
819  }
820  if (bx.smallEnd(1) == domain.smallEnd(1)) {
821  Box lo_y_dom_face(bx); lo_y_dom_face.setBig(1,bx.smallEnd(1));
822  if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir) {
823  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
824  rho_v_rhs(i,j,k) = zero;
825  });
826  } else if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir_upwind) {
827  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
828  if (v(i,j,k) >= zero) {
829  rho_v_rhs(i,j,k) = zero;
830  }
831  });
832  }
833  }
834  if (bx.bigEnd(1) == domain.bigEnd(1)) {
835  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);
836  if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir) {
837  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
838  rho_v_rhs(i,j,k) = zero;
839  });
840  } else if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir_upwind) {
841  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
842  if (v(i,j,k) <= zero) {
843  rho_v_rhs(i,j,k) = zero;
844  }
845  });
846  }
847  }
848 
849  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
850  { // z-momentum equation
851 
852  Real gpz = gpz_arr(i,j,k);
853 
854  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i,j,k-1)) : zero;
855 
856  rho_w_rhs(i, j, k) += (-gpz - abl_pressure_grad[2] + buoyancy_arr(i,j,k)) / (one + q) + zmom_src_arr(i,j,k);
857 
858  if (l_moving_terrain) {
859  rho_w_rhs(i, j, k) *= myhalf * (detJ_arr(i,j,k) + detJ_arr(i,j,k-1));
860  }
861  });
862 
863  auto const lo = lbound(bx);
864  auto const hi = ubound(bx);
865 
866  // Note: the logic below assumes no tiling in z!
867  if (level > 0) {
868 
869  const Array4<const Real>& rho_w_rhs_crse = zmom_crse_rhs->const_array(mfi);
870 
871  Box b2d = bx; b2d.setRange(2,0);
872 
873  if (lo.z > klo) {
874  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // bottom of box but not of domain
875  {
876  rho_w_rhs(i,j,lo.z) = rho_w_rhs_crse(i,j,lo.z);
877  });
878  }
879 
880  if (hi.z < khi+1) {
881  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // top of box but not of domain
882  {
883  rho_w_rhs(i,j,hi.z+1) = rho_w_rhs_crse(i,j,hi.z+1);
884  });
885  }
886  }
887 
888  {
889  BL_PROFILE("slow_rhs_pre_fluxreg");
890  // We only add to the flux registers in the final RK step
891  // NOTE: for now we are only refluxing density not (rho theta) since the latter seems to introduce
892  // a problem at top and bottom boundaries
893  if (l_reflux) {
894  int strt_comp_reflux = (l_fixed_rho) ? 1 : 0;
895  int num_comp_reflux = 1;
896  if (level < finest_level) {
897  fr_as_crse->CrseAdd(mfi,
898  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
899  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
900  }
901  if (level > 0) {
902  fr_as_fine->FineAdd(mfi,
903  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
904  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
905  }
906 
907  } // two-way coupling
908  } // end profile
909  } // mfi
910  } // OMP
911 }
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 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 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 > &u_tau_eb13, const amrex::Array4< const amrex::Real > &u_tau_eb23, const amrex::Array4< const amrex::Real > &v_tau_eb13, const amrex::Array4< const amrex::Real > &v_tau_eb23, const amrex::Array4< const amrex::Real > &w_tau_eb13, const amrex::Array4< const amrex::Real > &w_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 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
@ yface
Definition: ERF_EBStruct.H:20
@ zface
Definition: ERF_EBStruct.H:20
@ xface
Definition: ERF_EBStruct.H:20
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:255
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)
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
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
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:103
@ xvel_bc
Definition: ERF_IndexDefines.H:102
@ ext_dir
Definition: ERF_IndexDefines.H:243
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:251
@ gpz
Definition: ERF_IndexDefines.H:186
@ gpy
Definition: ERF_IndexDefines.H:185
@ gpx
Definition: ERF_IndexDefines.H:184
@ NumTypes
Definition: ERF_IndexDefines.H:196
@ ymom
Definition: ERF_IndexDefines.H:194
@ cons
Definition: ERF_IndexDefines.H:192
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193
@ qt
Definition: ERF_Kessler.H:28
@ ng
Definition: ERF_Morrison.H:48
@ xvel
Definition: ERF_IndexDefines.H:175
@ zvel
Definition: ERF_IndexDefines.H:177
@ yvel
Definition: ERF_IndexDefines.H:176
@ q
Definition: ERF_WSM6.H:169
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
Definition: ERF_EBStruct.H:27
EBBoundaryType eb_boundary_type
Definition: ERF_EBStruct.H:59
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1143
amrex::Real gravity
Definition: ERF_DataStruct.H:1220
bool use_shoc
Definition: ERF_DataStruct.H:1253
amrex::Vector< amrex::Real > vert_implicit_fac
Definition: ERF_DataStruct.H:1169
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > abl_pressure_grad
Definition: ERF_DataStruct.H:1306
bool implicit_thermal_diffusion
Definition: ERF_DataStruct.H:1171
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1146
amrex::Vector< int > anelastic
Definition: ERF_DataStruct.H:1152
AdvChoice advChoice
Definition: ERF_DataStruct.H:1142
MoistureType moisture_type
Definition: ERF_DataStruct.H:1299
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1125
amrex::Vector< int > fixed_density
Definition: ERF_DataStruct.H:1153
bool use_rotate_surface_flux
Definition: ERF_DataStruct.H:1250
EBChoice ebChoice
Definition: ERF_DataStruct.H:1147
CouplingType coupling_type
Definition: ERF_DataStruct.H:1298
Definition: ERF_TurbStruct.H:82
PBLType pbl_type
Definition: ERF_TurbStruct.H:459
bool use_keqn
Definition: ERF_TurbStruct.H:466
RANSType rans_type
Definition: ERF_TurbStruct.H:454
LESType les_type
Definition: ERF_TurbStruct.H:412
bool use_kturb
Definition: ERF_TurbStruct.H:465
Here is the call graph for this function: