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