ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_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, Vector< MultiFab > &S_scratch, 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, MultiFab *SmnSmn, MultiFab *eddyDiffs, MultiFab *Hfx1, MultiFab *Hfx2, MultiFab *Hfx3, MultiFab *Q1fx1, MultiFab *Q1fx2, MultiFab *Q1fx3, MultiFab *Q2fx3, MultiFab *Diss, 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 &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,
Vector< MultiFab > &  S_scratch,
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,
MultiFab *  SmnSmn,
MultiFab *  eddyDiffs,
MultiFab *  Hfx1,
MultiFab *  Hfx2,
MultiFab *  Hfx3,
MultiFab *  Q1fx1,
MultiFab *  Q1fx2,
MultiFab *  Q1fx3,
MultiFab *  Q2fx3,
MultiFab *  Diss,
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 &  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]S_scratchscratch space
[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
106 {
107  BL_PROFILE_REGION("erf_slow_rhs_pre()");
108 
109  const BCRec* bc_ptr_d = domain_bcs_type_d.data();
110  const BCRec* bc_ptr_h = domain_bcs_type_h.data();
111 
112  DiffChoice dc = solverChoice.diffChoice;
113  TurbChoice tc = solverChoice.turbChoice[level];
114 
115  const MultiFab* t_mean_mf = nullptr;
116  if (SurfLayer) { t_mean_mf = SurfLayer->get_mac_avg(level,2); }
117 
118  const Box& domain = geom.Domain();
119  int klo = domain.smallEnd(2);
120  int khi = domain.bigEnd(2);
121 
122  const AdvType l_horiz_adv_type = solverChoice.advChoice.dycore_horiz_adv_type;
123  const AdvType l_vert_adv_type = solverChoice.advChoice.dycore_vert_adv_type;
124  const Real l_horiz_upw_frac = solverChoice.advChoice.dycore_horiz_upw_frac;
125  const Real l_vert_upw_frac = solverChoice.advChoice.dycore_vert_upw_frac;
126  const bool l_use_terrain_fitted_coords = (solverChoice.mesh_type != MeshType::ConstantDz);
127  const bool l_moving_terrain = (solverChoice.terrain_type == TerrainType::MovingFittedMesh);
128  if (l_moving_terrain) AMREX_ALWAYS_ASSERT (l_use_terrain_fitted_coords);
129 
130  const bool l_use_mono_adv = solverChoice.use_mono_adv;
131  const bool l_reflux = ( (solverChoice.coupling_type == CouplingType::TwoWay) && (nrk == 2) && (finest_level > 0) );
132 
133  const bool l_use_diff = ( (dc.molec_diff_type != MolecDiffType::None) ||
134  (tc.les_type != LESType::None) ||
135  (tc.rans_type != RANSType::None) ||
136  (tc.pbl_type != PBLType::None) );
137  const bool l_use_turb = tc.use_kturb;
138  const bool l_need_SmnSmn = tc.use_keqn;
139 
140  const bool l_use_moisture = (solverChoice.moisture_type != MoistureType::None);
141  const bool l_use_SurfLayer = (SurfLayer != nullptr);
142  const bool l_rotate = (solverChoice.use_rotate_surface_flux);
143 
144  const bool l_anelastic = solverChoice.anelastic[level];
145  const bool l_fixed_rho = solverChoice.fixed_density;
146 
147  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
148  const Real* dx = geom.CellSize();
149 
150  // *****************************************************************************
151  // Combine external forcing terms
152  // *****************************************************************************
153  const Array<Real,AMREX_SPACEDIM> grav{0.0, 0.0, -solverChoice.gravity};
154  const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
155 
156  // *****************************************************************************
157  // Pre-computed quantities
158  // *****************************************************************************
159  int nvars = S_data[IntVars::cons].nComp();
160  const BoxArray& ba = S_data[IntVars::cons].boxArray();
161  const DistributionMapping& dm = S_data[IntVars::cons].DistributionMap();
162 
163  int nGhost = (solverChoice.terrain_type == TerrainType::EB) ? 2 : 1;
164  MultiFab Omega(convert(ba,IntVect(0,0,1)), dm, 1, nGhost);
165 
166  std::unique_ptr<MultiFab> expr;
167  std::unique_ptr<MultiFab> dflux_x;
168  std::unique_ptr<MultiFab> dflux_y;
169  std::unique_ptr<MultiFab> dflux_z;
170 
171  if (l_use_diff) {
172  erf_make_tau_terms(level,nrk,domain_bcs_type_h,z_phys_nd,
173  S_data,xvel,yvel,zvel,Tau_lev,
174  SmnSmn,eddyDiffs,geom,solverChoice,SurfLayer,
175  detJ,mapfac);
176 
177  dflux_x = std::make_unique<MultiFab>(convert(ba,IntVect(1,0,0)), dm, nvars, 0);
178  dflux_y = std::make_unique<MultiFab>(convert(ba,IntVect(0,1,0)), dm, nvars, 0);
179  dflux_z = std::make_unique<MultiFab>(convert(ba,IntVect(0,0,1)), dm, nvars, 0);
180 
181  if (l_use_SurfLayer) {
182  Vector<const MultiFab*> mfs = {&S_data[IntVars::cons], &xvel, &yvel, &zvel};
183  SurfLayer->impose_SurfaceLayer_bcs(level, mfs, Tau_lev,
184  Hfx1, Hfx2, Hfx3,
185  Q1fx1, Q1fx2, Q1fx3,
186  &z_phys_nd);
187  }
188  } // l_use_diff
189 
190  // *****************************************************************************
191  // Monotonic advection for scalars
192  // *****************************************************************************
193  int nvar = S_data[IntVars::cons].nComp();
194  Vector<Real> max_scal(nvar, 1.0e34); Gpu::DeviceVector<Real> max_scal_d(nvar);
195  Vector<Real> min_scal(nvar,-1.0e34); Gpu::DeviceVector<Real> min_scal_d(nvar);
196  if (l_use_mono_adv) {
197  auto const& ma_s_arr = S_data[IntVars::cons].const_arrays();
198  for (int ivar(RhoTheta_comp); ivar<RhoKE_comp; ++ivar) {
199  GpuTuple<Real,Real> mm = ParReduce(TypeList<ReduceOpMax,ReduceOpMin>{},
200  TypeList<Real, Real>{},
201  S_data[IntVars::cons], IntVect(0),
202  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
203  -> GpuTuple<Real,Real>
204  {
205  return { ma_s_arr[box_no](i,j,k,ivar), ma_s_arr[box_no](i,j,k,ivar) };
206  });
207  max_scal[ivar] = get<0>(mm);
208  min_scal[ivar] = get<1>(mm);
209  }
210  }
211  Gpu::copy(Gpu::hostToDevice, max_scal.begin(), max_scal.end(), max_scal_d.begin());
212  Gpu::copy(Gpu::hostToDevice, min_scal.begin(), min_scal.end(), min_scal_d.begin());
213  Real* max_s_ptr = max_scal_d.data();
214  Real* min_s_ptr = min_scal_d.data();
215 
216  // This is just cautionary to deal with grid boundaries that aren't domain boundaries
217  S_rhs[IntVars::zmom].setVal(0.0);
218 
219  // *****************************************************************************
220  // Define updates and fluxes in the current RK stage
221  // *****************************************************************************
222 #ifdef _OPENMP
223 #pragma omp parallel if (Gpu::notInLaunchRegion())
224 #endif
225  {
226  std::array<FArrayBox,AMREX_SPACEDIM> flux;
227  std::array<FArrayBox,AMREX_SPACEDIM> flux_u;
228  std::array<FArrayBox,AMREX_SPACEDIM> flux_v;
229  std::array<FArrayBox,AMREX_SPACEDIM> flux_w;
230  std::array<FArrayBox,AMREX_SPACEDIM> flux_tmp;
231 
232  // Cell-centered masks for EB (used for flux interpolation)
233  bool already_on_centroids = false;
234  Vector<iMultiFab> physbnd_mask;
235  physbnd_mask.resize(IntVars::NumTypes);
236  if (solverChoice.terrain_type == TerrainType::EB) {
237  physbnd_mask[IntVars::cons].define(S_data[IntVars::cons].boxArray(), S_data[IntVars::cons].DistributionMap(), 1, 1);
238  physbnd_mask[IntVars::cons].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
239  // physbnd_mask[IntVars::cons].FillBoundary(geom.periodicity());
240  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
241  physbnd_mask[1+dir].define(S_data[1+dir].boxArray(), S_data[1+dir].DistributionMap(), 1, 1);
242  physbnd_mask[1+dir].BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
243  // physbnd_mask[1+dir].FillBoundary(geom.periodicity());
244  }
245  }
246 
247  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
248  {
249  Box bx = mfi.tilebox();
250  Box tbx = mfi.nodaltilebox(0);
251  Box tby = mfi.nodaltilebox(1);
252  Box tbz = mfi.nodaltilebox(2);
253 
254  // Boxes for momentum fluxes
255  Vector<Box> tbx_grown(AMREX_SPACEDIM);
256  Vector<Box> tby_grown(AMREX_SPACEDIM);
257  Vector<Box> tbz_grown(AMREX_SPACEDIM);
258  if (solverChoice.terrain_type == TerrainType::EB) {
259  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
260  tbx_grown[dir] = tbx;
261  tby_grown[dir] = tby;
262  tbz_grown[dir] = tbz;
263  IntVect iv(1, 1, 1);
264  iv[dir] = 0;
265  tbx_grown[dir] = (tbx_grown[dir].growHi(dir,1)).grow(iv);
266  tby_grown[dir] = (tby_grown[dir].growHi(dir,1)).grow(iv);
267  tbz_grown[dir] = (tbz_grown[dir].growHi(dir,1)).grow(iv);
268  }
269  }
270 
271  // We don't compute a source term for z-momentum on the bottom or top domain boundary
272  if (tbz.smallEnd(2) == domain.smallEnd(2)) {
273  tbz.growLo(2,-1);
274  }
275  if (tbz.bigEnd(2) == domain.bigEnd(2)+1) {
276  tbz.growHi(2,-1);
277  }
278 
279  const Array4<const Real> & cell_data = S_data[IntVars::cons].array(mfi);
280  const Array4<const Real> & cell_prim = S_prim.array(mfi);
281  const Array4<Real> & cell_rhs = S_rhs[IntVars::cons].array(mfi);
282 
283  const Array4<const Real> & cell_old = S_old[IntVars::cons].array(mfi);
284 
285  const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
286  const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
287  const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
288  const Array4<Real const>& buoyancy_arr = buoyancy.const_array(mfi);
289 
290  const Array4<Real const>& gpx_arr = gradp[GpVars::gpx].const_array(mfi);
291  const Array4<Real const>& gpy_arr = gradp[GpVars::gpy].const_array(mfi);
292  const Array4<Real const>& gpz_arr = gradp[GpVars::gpz].const_array(mfi);
293 
294  const Array4<Real const>& qt_arr = qt.const_array(mfi);
295 
296  const Array4<Real>& rho_u_old = S_old[IntVars::xmom].array(mfi);
297  const Array4<Real>& rho_v_old = S_old[IntVars::ymom].array(mfi);
298 
299  if (l_anelastic) {
300  // When anelastic we must reset these to 0 each RK step
301  S_scratch[IntVars::xmom][mfi].template setVal<RunOn::Device>(0.0,tbx);
302  S_scratch[IntVars::ymom][mfi].template setVal<RunOn::Device>(0.0,tby);
303  S_scratch[IntVars::zmom][mfi].template setVal<RunOn::Device>(0.0,tbz);
304  }
305 
306  Array4<Real> avg_xmom = S_scratch[IntVars::xmom].array(mfi);
307  Array4<Real> avg_ymom = S_scratch[IntVars::ymom].array(mfi);
308  Array4<Real> avg_zmom = S_scratch[IntVars::zmom].array(mfi);
309 
310  const Array4<const Real> & u = xvel.array(mfi);
311  const Array4<const Real> & v = yvel.array(mfi);
312  const Array4<const Real> & w = zvel.array(mfi);
313 
314  const Array4<const Real>& rho_u = S_data[IntVars::xmom].array(mfi);
315  const Array4<const Real>& rho_v = S_data[IntVars::ymom].array(mfi);
316  const Array4<const Real>& rho_w = S_data[IntVars::zmom].array(mfi);
317 
318  // Map factors
319  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
320  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
321  const Array4<const Real>& mf_vx = mapfac[MapFacType::v_x]->const_array(mfi);
322  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
323  const Array4<const Real>& mf_uy = mapfac[MapFacType::u_y]->const_array(mfi);
324  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
325 
326  const Array4< Real>& omega_arr = Omega.array(mfi);
327 
328  Array4<const Real> z_t;
329  if (z_t_mf) {
330  z_t = z_t_mf->array(mfi);
331  } else {
332  z_t = Array4<const Real>{};
333  }
334 
335  const Array4<Real>& rho_u_rhs = S_rhs[IntVars::xmom].array(mfi);
336  const Array4<Real>& rho_v_rhs = S_rhs[IntVars::ymom].array(mfi);
337  const Array4<Real>& rho_w_rhs = S_rhs[IntVars::zmom].array(mfi);
338 
339  const Array4<Real const>& mu_turb = l_use_turb ? eddyDiffs->const_array(mfi) : Array4<const Real>{};
340 
341  // Terrain metrics
342  const Array4<const Real>& z_nd = z_phys_nd.const_array(mfi);
343 
344  // *****************************************************************************
345  // Define flux arrays for use in advection
346  // *****************************************************************************
347  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
348  if (solverChoice.terrain_type != TerrainType::EB) {
349  flux[dir].resize(surroundingNodes(bx,dir),2);
350  } else {
351  flux[dir].resize(surroundingNodes(bx,dir).grow(1),2);
352  }
353  flux[dir].setVal<RunOn::Device>(0.);
354  if (l_use_mono_adv) {
355  flux_tmp[dir].resize(surroundingNodes(bx,dir),2);
356  flux_tmp[dir].setVal<RunOn::Device>(0.);
357  }
358  }
359  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
360  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
361 
362  // Define flux arrays for momentum variables (used only for EB now)
363  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_u_arr{};
364  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_v_arr{};
365  GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_w_arr{};
366  if (solverChoice.terrain_type == TerrainType::EB) {
367  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
368  flux_u[dir].resize(tbx_grown[dir],1);
369  flux_v[dir].resize(tby_grown[dir],1);
370  flux_w[dir].resize(tbz_grown[dir],1);
371  flux_u[dir].setVal<RunOn::Device>(0.);
372  flux_v[dir].setVal<RunOn::Device>(0.);
373  flux_w[dir].setVal<RunOn::Device>(0.);
374  flx_u_arr[dir] = flux_u[dir].array();
375  flx_v_arr[dir] = flux_v[dir].array();
376  flx_w_arr[dir] = flux_w[dir].array();
377  }
378  }
379 
380  Array4<Real> tmpx = (l_use_mono_adv) ? flux_tmp[0].array() : Array4<Real>{};
381  Array4<Real> tmpy = (l_use_mono_adv) ? flux_tmp[1].array() : Array4<Real>{};
382  Array4<Real> tmpz = (l_use_mono_adv) ? flux_tmp[2].array() : Array4<Real>{};
383  const GpuArray<Array4<Real>, AMREX_SPACEDIM> flx_tmp_arr{{AMREX_D_DECL(tmpx,tmpy,tmpz)}};
384 
385  // *****************************************************************************
386  // Contravariant flux field
387  // *****************************************************************************
388  {
389  BL_PROFILE("slow_rhs_making_omega");
390  IntVect nGrowVect = (solverChoice.terrain_type == TerrainType::EB)
391  ? IntVect(AMREX_D_DECL(2, 2, 2)) : IntVect(AMREX_D_DECL(1, 1, 1));
392  Box gbxo = surroundingNodes(bx,2); gbxo.grow(nGrowVect);
393  //
394  // Now create Omega with momentum (not velocity) with z_t subtracted if moving terrain
395  // ONLY if not doing anelastic + terrain -- in that case Omega will be defined coming
396  // out of the projection
397  //
398  if (!l_use_terrain_fitted_coords) {
399  ParallelFor(gbxo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
400  omega_arr(i,j,k) = rho_w(i,j,k);
401  });
402 
403  } else {
404 
405  Box gbxo_lo = gbxo; gbxo_lo.setBig(2,domain.smallEnd(2));
406  int lo_z_face = domain.smallEnd(2);
407  if (gbxo_lo.smallEnd(2) <= lo_z_face) {
408  ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
409  omega_arr(i,j,k) = 0.;
410  });
411  }
412  Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
413  int hi_z_face = domain.bigEnd(2)+1;
414  if (gbxo_hi.bigEnd(2) >= hi_z_face) {
415  ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
416  omega_arr(i,j,k) = rho_w(i,j,k);
417  });
418  }
419 
420  if (z_t) { // Note we never do anelastic with moving terrain
421  Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
422  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
423  // We define rho on the z-face the same way as in MomentumToVelocity/VelocityToMomentum
424  Real rho_at_face = 0.5 * (cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp));
425  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
426  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv) -
427  rho_at_face * z_t(i,j,k);
428  });
429  } else {
430  Box gbxo_mid = gbxo;
431  if (gbxo_mid.smallEnd(2) <= domain.smallEnd(2)) {
432  gbxo_mid.setSmall(2,1);
433  }
434  if (gbxo_mid.bigEnd(2) >= domain.bigEnd(2)+1) {
435  gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
436  }
437  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
438  omega_arr(i,j,k) = OmegaFromW(i,j,k,rho_w(i,j,k),
439  rho_u,rho_v,mf_ux,mf_vy,z_nd,dxInv);
440  });
441  }
442  }
443  } // end profile
444 
445  // *****************************************************************************
446  // Diffusive terms (pre-computed above)
447  // *****************************************************************************
448  // No terrain diffusion
449  Array4<Real> tau11,tau22,tau33;
450  Array4<Real> tau12,tau13,tau23;
451  if (Tau_lev[TauType::tau11]) {
452  tau11 = Tau_lev[TauType::tau11]->array(mfi); tau22 = Tau_lev[TauType::tau22]->array(mfi);
453  tau33 = Tau_lev[TauType::tau33]->array(mfi); tau12 = Tau_lev[TauType::tau12]->array(mfi);
454  tau13 = Tau_lev[TauType::tau13]->array(mfi); tau23 = Tau_lev[TauType::tau23]->array(mfi);
455  } else {
456  tau11 = Array4<Real>{}; tau22 = Array4<Real>{}; tau33 = Array4<Real>{};
457  tau12 = Array4<Real>{}; tau13 = Array4<Real>{}; tau23 = Array4<Real>{};
458  }
459  // Terrain diffusion
460  Array4<Real> tau21,tau31,tau32;
461  if (Tau_lev[TauType::tau21]) {
462  tau21 = Tau_lev[TauType::tau21]->array(mfi);
463  tau31 = Tau_lev[TauType::tau31]->array(mfi);
464  tau32 = Tau_lev[TauType::tau32]->array(mfi);
465  } else {
466  tau21 = Array4<Real>{}; tau31 = Array4<Real>{}; tau32 = Array4<Real>{};
467  }
468 
469  // Strain magnitude
470  Array4<Real> SmnSmn_a;
471  if (l_need_SmnSmn) {
472  SmnSmn_a = SmnSmn->array(mfi);
473  } else {
474  SmnSmn_a = Array4<Real>{};
475  }
476 
477  // *****************************************************************************
478  // Define updates in the RHS of continuity and potential temperature equations
479  // *****************************************************************************
480  Array4<const int> mask_arr{};
481  Array4<const EBCellFlag> cfg_arr{};
482  Array4<const Real> ax_arr{};
483  Array4<const Real> ay_arr{};
484  Array4<const Real> az_arr{};
485  Array4<const Real> fcx_arr{};
486  Array4<const Real> fcy_arr{};
487  Array4<const Real> fcz_arr{};
488  Array4<const Real> detJ_arr{};
489  if (solverChoice.terrain_type == TerrainType::EB)
490  {
491  EBCellFlagFab const& cfg = (ebfact.get_const_factory())->getMultiEBCellFlagFab()[mfi];
492  cfg_arr = cfg.const_array();
493  ax_arr = (ebfact.get_const_factory())->getAreaFrac()[0]->const_array(mfi);
494  ay_arr = (ebfact.get_const_factory())->getAreaFrac()[1]->const_array(mfi);
495  az_arr = (ebfact.get_const_factory())->getAreaFrac()[2]->const_array(mfi);
496  fcx_arr = (ebfact.get_const_factory())->getFaceCent()[0]->const_array(mfi);
497  fcy_arr = (ebfact.get_const_factory())->getFaceCent()[1]->const_array(mfi);
498  fcz_arr = (ebfact.get_const_factory())->getFaceCent()[2]->const_array(mfi);
499  detJ_arr = (ebfact.get_const_factory())->getVolFrac().const_array(mfi);
500  // if (!already_on_centroids) {mask_arr = physbnd_mask[IntVars::cons].const_array(mfi);}
501  mask_arr = physbnd_mask[IntVars::cons].const_array(mfi);
502  } else {
503  ax_arr = ax.const_array(mfi);
504  ay_arr = ay.const_array(mfi);
505  az_arr = az.const_array(mfi);
506  detJ_arr = detJ.const_array(mfi);
507  }
508 
509  if (solverChoice.terrain_type != TerrainType::EB){
510  AdvectionSrcForRho(bx, cell_rhs,
511  rho_u, rho_v, omega_arr, // these are being used to build the fluxes
512  avg_xmom, avg_ymom, avg_zmom, // these are being defined from the fluxes
513  ax_arr, ay_arr, az_arr, detJ_arr,
514  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
515  flx_arr, l_fixed_rho);
516  } else {
517  EBAdvectionSrcForRho(bx, cell_rhs,
518  rho_u, rho_v, omega_arr,
519  avg_xmom, avg_ymom, avg_zmom,
520  mask_arr, cfg_arr,
521  ax_arr, ay_arr, az_arr,
522  fcx_arr, fcy_arr, fcz_arr, detJ_arr,
523  dxInv, mf_mx, mf_my, mf_uy, mf_vx,
524  flx_arr, l_fixed_rho,
525  already_on_centroids);
526  }
527 
528  int icomp = RhoTheta_comp; int ncomp = 1;
529  if (solverChoice.terrain_type != TerrainType::EB){
530  AdvectionSrcForScalars(dt, bx, icomp, ncomp,
531  avg_xmom, avg_ymom, avg_zmom,
532  cell_data, cell_prim, cell_rhs,
533  l_use_mono_adv, max_s_ptr, min_s_ptr,
534  detJ_arr, dxInv, mf_mx, mf_my,
535  l_horiz_adv_type, l_vert_adv_type,
536  l_horiz_upw_frac, l_vert_upw_frac,
537  flx_arr, flx_tmp_arr, domain, bc_ptr_h);
538  } else {
539  EBAdvectionSrcForScalars(bx, icomp, ncomp,
540  avg_xmom, avg_ymom, avg_zmom,
541  cell_prim, cell_rhs,
542  mask_arr, cfg_arr, ax_arr, ay_arr, az_arr,
543  fcx_arr, fcy_arr, fcz_arr,
544  detJ_arr, dxInv, mf_mx, mf_my,
545  l_horiz_adv_type, l_vert_adv_type,
546  l_horiz_upw_frac, l_vert_upw_frac,
547  flx_arr, domain, bc_ptr_h,
548  already_on_centroids);
549  }
550 
551  if (l_use_diff) {
552  Array4<Real> diffflux_x = dflux_x->array(mfi);
553  Array4<Real> diffflux_y = dflux_y->array(mfi);
554  Array4<Real> diffflux_z = dflux_z->array(mfi);
555 
556  Array4<Real> hfx_x = Hfx1->array(mfi);
557  Array4<Real> hfx_y = Hfx2->array(mfi);
558  Array4<Real> hfx_z = Hfx3->array(mfi);
559 
560  Array4<Real> q1fx_x = (Q1fx1) ? Q1fx1->array(mfi) : Array4<Real>{};
561  Array4<Real> q1fx_y = (Q1fx2) ? Q1fx2->array(mfi) : Array4<Real>{};
562  Array4<Real> q1fx_z = (Q1fx3) ? Q1fx3->array(mfi) : Array4<Real>{};
563 
564  Array4<Real> q2fx_z = (Q2fx3) ? Q2fx3->array(mfi) : Array4<Real>{};
565  Array4<Real> diss = Diss->array(mfi);
566 
567  const Array4<const Real> tm_arr = t_mean_mf ? t_mean_mf->const_array(mfi) : Array4<const Real>{};
568 
569  // NOTE: No diffusion for continuity, so n starts at 1.
570  int n_start = RhoTheta_comp;
571  int n_comp = 1;
572 
573  if (l_use_terrain_fitted_coords) {
574  DiffusionSrcForState_T(bx, domain, n_start, n_comp, l_rotate, u, v,
575  cell_data, cell_prim, cell_rhs,
576  diffflux_x, diffflux_y, diffflux_z,
577  z_nd, ax_arr, ay_arr, az_arr, detJ_arr,
578  dxInv, SmnSmn_a,
579  mf_mx, mf_ux, mf_vx,
580  mf_my, mf_uy, mf_vy,
581  hfx_x, hfx_y, hfx_z, q1fx_x, q1fx_y, q1fx_z, q2fx_z, diss,
582  mu_turb, solverChoice, level,
583  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer);
584  } else {
585  DiffusionSrcForState_N(bx, domain, n_start, n_comp, u, v,
586  cell_data, cell_prim, cell_rhs,
587  diffflux_x, diffflux_y, diffflux_z,
588  dxInv, SmnSmn_a,
589  mf_mx, mf_ux, mf_vx,
590  mf_my, mf_uy, mf_vy,
591  hfx_z, q1fx_z, q2fx_z, diss,
592  mu_turb, solverChoice, level,
593  tm_arr, grav_gpu, bc_ptr_d, l_use_SurfLayer);
594  }
595  }
596 
597  const Array4<Real const>& source_arr = cc_src.const_array(mfi);
598  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
599  {
600  cell_rhs(i,j,k,Rho_comp) += source_arr(i,j,k,Rho_comp);
601  cell_rhs(i,j,k,RhoTheta_comp) += source_arr(i,j,k,RhoTheta_comp);
602  });
603 
604  // Multiply the slow RHS for rho and rhotheta by detJ here so we don't have to later
605  if (l_use_terrain_fitted_coords && l_moving_terrain) {
606  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
607  {
608  cell_rhs(i,j,k,Rho_comp) *= detJ_arr(i,j,k);
609  cell_rhs(i,j,k,RhoTheta_comp) *= detJ_arr(i,j,k);
610  });
611  }
612 
613  // If anelastic and in second RK stage, take average of old-time and new-time source
614  if ( l_anelastic && (nrk == 1) )
615  {
616  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
617  {
618  cell_rhs(i,j,k, Rho_comp) *= 0.5;
619  cell_rhs(i,j,k,RhoTheta_comp) *= 0.5;
620 
621  cell_rhs(i,j,k, Rho_comp) += 0.5 / dt * (cell_data(i,j,k, Rho_comp) - cell_old(i,j,k, Rho_comp));
622  cell_rhs(i,j,k,RhoTheta_comp) += 0.5 / dt * (cell_data(i,j,k,RhoTheta_comp) - cell_old(i,j,k,RhoTheta_comp));
623  });
624  }
625 
626  // *****************************************************************************
627  // Define updates in the RHS of {x, y, z}-momentum equations
628  // *****************************************************************************
629  int lo_z_face = domain.smallEnd(2);
630  int hi_z_face = domain.bigEnd(2)+1;
631 
632  AdvectionSrcForMom(mfi, bx, tbx, tby, tbz, tbx_grown, tby_grown, tbz_grown,
633  rho_u_rhs, rho_v_rhs, rho_w_rhs,
634  cell_data, u, v, w,
635  rho_u, rho_v, omega_arr,
636  z_nd, ax_arr, ay_arr, az_arr, detJ_arr,
637  stretched_dz_d,
638  dxInv, mf_mx, mf_ux, mf_vx, mf_my, mf_uy, mf_vy,
639  l_horiz_adv_type, l_vert_adv_type,
640  l_horiz_upw_frac, l_vert_upw_frac,
641  solverChoice.mesh_type, solverChoice.terrain_type,
642  ebfact, flx_u_arr, flx_v_arr, flx_w_arr,
643  physbnd_mask, already_on_centroids,
644  lo_z_face, hi_z_face, domain, bc_ptr_h);
645 
646  if (l_use_diff) {
647  // Note: tau** were calculated with calls to
648  // ComputeStress[Cons|Var]Visc_[N|T] in which ConsVisc ("constant
649  // viscosity") means that there is no contribution from a
650  // turbulence model. However, whether this field truly is constant
651  // depends on whether MolecDiffType is Constant or ConstantAlpha.
652  if (l_use_terrain_fitted_coords) {
653  DiffusionSrcForMom_T(tbx, tby, tbz,
654  rho_u_rhs, rho_v_rhs, rho_w_rhs,
655  tau11, tau22, tau33,
656  tau12, tau21,
657  tau13, tau31,
658  tau23, tau32,
659  detJ_arr, dxInv,
660  mf_mx, mf_ux, mf_vx,
661  mf_my, mf_uy, mf_vy);
662  } else {
663  DiffusionSrcForMom_N(tbx, tby, tbz,
664  rho_u_rhs, rho_v_rhs, rho_w_rhs,
665  tau11, tau22, tau33,
666  tau12, tau13, tau23, dxInv,
667  mf_mx, mf_ux, mf_vx,
668  mf_my, mf_uy, mf_vy);
669  }
670  }
671 
672  auto abl_pressure_grad = solverChoice.abl_pressure_grad;
673 
674  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
675  { // x-momentum equation
676 
677  Real gpx = gpx_arr(i,j,k) * mf_ux(i,j,0);
678 
679  Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) : 0.0;
680 
681  rho_u_rhs(i, j, k) += (-gpx - abl_pressure_grad[0]) / (1.0 + q) + xmom_src_arr(i,j,k);
682 
683  if (l_moving_terrain) {
684  Real h_zeta = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd);
685  rho_u_rhs(i, j, k) *= h_zeta;
686  }
687 
688  if ( l_anelastic && (nrk == 1) ) {
689  rho_u_rhs(i,j,k) *= 0.5;
690  rho_u_rhs(i,j,k) += 0.5 / dt * (rho_u(i,j,k) - rho_u_old(i,j,k));
691  }
692  });
693 
694  ParallelFor(tby, [=] AMREX_GPU_DEVICE (int i, int j, int k)
695  { // y-momentum equation
696 
697  Real gpy = gpy_arr(i,j,k) * mf_vy(i,j,0);
698 
699  Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) : 0.0;
700 
701  rho_v_rhs(i, j, k) += (-gpy - abl_pressure_grad[1]) / (1.0 + q) + ymom_src_arr(i,j,k);
702 
703  if (l_moving_terrain) {
704  Real h_zeta = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd);
705  rho_v_rhs(i, j, k) *= h_zeta;
706  }
707 
708  if ( l_anelastic && (nrk == 1) ) {
709  rho_v_rhs(i,j,k) *= 0.5;
710  rho_v_rhs(i,j,k) += 0.5 / dt * (rho_v(i,j,k) - rho_v_old(i,j,k));
711  }
712  });
713 
714  // *****************************************************************************
715  // Zero out source terms for x- and y- momenta if at walls or inflow
716  // We need to do this -- even though we call the boundary conditions later --
717  // because the slow source is used to update the state in the fast interpolater.
718  // *****************************************************************************
719  if (bx.smallEnd(0) == domain.smallEnd(0)) {
720  Box lo_x_dom_face(bx); lo_x_dom_face.setBig(0,bx.smallEnd(0));
721  if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir) {
722  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
723  rho_u_rhs(i,j,k) = 0.;
724  });
725  } else if (bc_ptr_h[BCVars::xvel_bc].lo(0) == ERFBCType::ext_dir_upwind) {
726  ParallelFor(lo_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
727  if (u(i,j,k) >= 0.) {
728  rho_u_rhs(i,j,k) = 0.;
729  }
730  });
731  }
732  }
733  if (bx.bigEnd(0) == domain.bigEnd(0)) {
734  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);
735  if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir) {
736  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
737  rho_u_rhs(i,j,k) = 0.;
738  });
739  } else if (bc_ptr_h[BCVars::xvel_bc].hi(0) == ERFBCType::ext_dir_upwind) {
740  ParallelFor(hi_x_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
741  if (u(i,j,k) <= 0.) {
742  rho_u_rhs(i,j,k) = 0.;
743  }
744  });
745  }
746  }
747  if (bx.smallEnd(1) == domain.smallEnd(1)) {
748  Box lo_y_dom_face(bx); lo_y_dom_face.setBig(1,bx.smallEnd(1));
749  if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir) {
750  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
751  rho_v_rhs(i,j,k) = 0.;
752  });
753  } else if (bc_ptr_h[BCVars::yvel_bc].lo(1) == ERFBCType::ext_dir_upwind) {
754  ParallelFor(lo_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
755  if (v(i,j,k) >= 0.) {
756  rho_v_rhs(i,j,k) = 0.;
757  }
758  });
759  }
760  }
761  if (bx.bigEnd(1) == domain.bigEnd(1)) {
762  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);
763  if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir) {
764  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
765  rho_v_rhs(i,j,k) = 0.;
766  });
767  } else if (bc_ptr_h[BCVars::yvel_bc].hi(1) == ERFBCType::ext_dir_upwind) {
768  ParallelFor(hi_y_dom_face, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
769  if (v(i,j,k) <= 0.) {
770  rho_v_rhs(i,j,k) = 0.;
771  }
772  });
773  }
774  }
775 
776  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
777  { // z-momentum equation
778 
779  Real gpz = gpz_arr(i,j,k);
780 
781  Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i,j,k-1)) : 0.0;
782 
783  rho_w_rhs(i, j, k) += (-gpz - abl_pressure_grad[2] + buoyancy_arr(i,j,k)) / (1.0 + q) + zmom_src_arr(i,j,k);
784 
785  if (l_use_terrain_fitted_coords && l_moving_terrain) {
786  rho_w_rhs(i, j, k) *= 0.5 * (detJ_arr(i,j,k) + detJ_arr(i,j,k-1));
787  }
788  });
789 
790  auto const lo = lbound(bx);
791  auto const hi = ubound(bx);
792 
793  // Note: the logic below assumes no tiling in z!
794  if (level > 0) {
795 
796  const Array4<const Real>& rho_w_rhs_crse = zmom_crse_rhs->const_array(mfi);
797 
798  Box b2d = bx; b2d.setRange(2,0);
799 
800  if (lo.z > klo) {
801  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // bottom of box but not of domain
802  {
803  rho_w_rhs(i,j,lo.z) = rho_w_rhs_crse(i,j,lo.z);
804  });
805  }
806 
807  if (hi.z < khi+1) {
808  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int ) // top of box but not of domain
809  {
810  rho_w_rhs(i,j,hi.z+1) = rho_w_rhs_crse(i,j,hi.z+1);
811  });
812  }
813  }
814 
815  {
816  BL_PROFILE("slow_rhs_pre_fluxreg");
817  // We only add to the flux registers in the final RK step
818  // NOTE: for now we are only refluxing density not (rho theta) since the latter seems to introduce
819  // a problem at top and bottom boundaries
820  if (l_reflux) {
821  int strt_comp_reflux = (l_fixed_rho) ? 1 : 0;
822  int num_comp_reflux = 1;
823  if (level < finest_level) {
824  fr_as_crse->CrseAdd(mfi,
825  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
826  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
827  }
828  if (level > 0) {
829  fr_as_fine->FineAdd(mfi,
830  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
831  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
832  }
833 
834  // This is necessary here so we don't go on to the next FArrayBox without
835  // having finished copying the fluxes into the FluxRegisters (since the fluxes
836  // are stored in temporary FArrayBox's)
837  Gpu::streamSynchronize();
838 
839  } // two-way coupling
840  } // end profile
841  } // mfi
842  } // OMP
843 }
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 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)
void AdvectionSrcForScalars(const amrex::Real &dt, 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 > &cur_cons, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< amrex::Real > &src, const bool &use_mono_adv, amrex::Real *max_s_ptr, amrex::Real *min_s_ptr, 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::GpuArray< amrex::Array4< amrex::Real >, AMREX_SPACEDIM > &flx_tmp_arr, const amrex::Box &domain, const amrex::BCRec *bc_ptr_h)
@ tau12
Definition: ERF_DataStruct.H:30
@ tau23
Definition: ERF_DataStruct.H:30
@ tau33
Definition: ERF_DataStruct.H:30
@ tau22
Definition: ERF_DataStruct.H:30
@ tau11
Definition: ERF_DataStruct.H:30
@ tau32
Definition: ERF_DataStruct.H:30
@ tau31
Definition: ERF_DataStruct.H:30
@ tau21
Definition: ERF_DataStruct.H:30
@ tau13
Definition: ERF_DataStruct.H:30
@ nvars
Definition: ERF_DataStruct.H:87
@ v_x
Definition: ERF_DataStruct.H:22
@ u_y
Definition: ERF_DataStruct.H:23
@ v_y
Definition: ERF_DataStruct.H:23
@ m_y
Definition: ERF_DataStruct.H:23
@ u_x
Definition: ERF_DataStruct.H:22
@ m_x
Definition: ERF_DataStruct.H:22
void DiffusionSrcForMom_N(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 > &tau13, const amrex::Array4< const amrex::Real > &tau23, 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)
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)
void DiffusionSrcForMom_T(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, 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)
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 > &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)
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)
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
AdvType
Definition: ERF_IndexDefines.H:221
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
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, MultiFab *SmnSmn, MultiFab *eddyDiffs, const Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &, const MultiFab &detJ, Vector< std::unique_ptr< MultiFab >> &mapfac)
Definition: ERF_MakeTauTerms.cpp:12
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:415
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:96
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:139
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:40
@ yvel_bc
Definition: ERF_IndexDefines.H:88
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ ext_dir
Definition: ERF_IndexDefines.H:209
@ ext_dir_upwind
Definition: ERF_IndexDefines.H:217
@ gpz
Definition: ERF_IndexDefines.H:152
@ gpy
Definition: ERF_IndexDefines.H:151
@ gpx
Definition: ERF_IndexDefines.H:150
@ NumTypes
Definition: ERF_IndexDefines.H:162
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
@ qt
Definition: ERF_Kessler.H:27
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
AdvType dycore_vert_adv_type
Definition: ERF_AdvStruct.H:340
amrex::Real dycore_vert_upw_frac
Definition: ERF_AdvStruct.H:350
AdvType dycore_horiz_adv_type
Definition: ERF_AdvStruct.H:339
amrex::Real dycore_horiz_upw_frac
Definition: ERF_AdvStruct.H:349
Definition: ERF_DiffStruct.H:19
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:81
static MeshType mesh_type
Definition: ERF_DataStruct.H:691
bool use_mono_adv
Definition: ERF_DataStruct.H:781
DiffChoice diffChoice
Definition: ERF_DataStruct.H:700
amrex::Real gravity
Definition: ERF_DataStruct.H:740
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > abl_pressure_grad
Definition: ERF_DataStruct.H:794
bool fixed_density
Definition: ERF_DataStruct.H:709
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:702
amrex::Vector< int > anelastic
Definition: ERF_DataStruct.H:707
AdvChoice advChoice
Definition: ERF_DataStruct.H:699
MoistureType moisture_type
Definition: ERF_DataStruct.H:787
static TerrainType terrain_type
Definition: ERF_DataStruct.H:685
bool use_rotate_surface_flux
Definition: ERF_DataStruct.H:768
CouplingType coupling_type
Definition: ERF_DataStruct.H:786
Definition: ERF_TurbStruct.H:31
PBLType pbl_type
Definition: ERF_TurbStruct.H:276
bool use_keqn
Definition: ERF_TurbStruct.H:283
RANSType rans_type
Definition: ERF_TurbStruct.H:273
LESType les_type
Definition: ERF_TurbStruct.H:236
bool use_kturb
Definition: ERF_TurbStruct.H:282
Here is the call graph for this function: