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