ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SlowRhsPost.cpp File Reference
#include <AMReX.H>
#include <ERF_SrcHeaders.H>
#include <ERF_TI_slow_headers.H>
#include <ERF_EBAdvection.H>
#include <ERF_EBRedistribute.H>
Include dependency graph for ERF_SlowRhsPost.cpp:

Functions

void erf_slow_rhs_post (int level, int finest_level, int nrk, Real dt, int n_qstate, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old, Vector< MultiFab > &S_new, Vector< MultiFab > &S_data, const MultiFab &S_prim, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &, const MultiFab &source, const MultiFab *SmnSmn, const 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, std::unique_ptr< MultiFab > &z_phys_nd, std::unique_ptr< MultiFab > &z_phys_cc, std::unique_ptr< MultiFab > &ax, std::unique_ptr< MultiFab > &ay, std::unique_ptr< MultiFab > &az, std::unique_ptr< MultiFab > &detJ, MultiFab *detJ_new, Gpu::DeviceVector< Real > &stretched_dz_d, Vector< std::unique_ptr< MultiFab >> &mapfac, amrex::EBFArrayBoxFactory const &ebfact, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine)
 

Function Documentation

◆ erf_slow_rhs_post()

void erf_slow_rhs_post ( int  level,
int  finest_level,
int  nrk,
Real  dt,
int  n_qstate,
Vector< MultiFab > &  S_rhs,
Vector< MultiFab > &  S_old,
Vector< MultiFab > &  S_new,
Vector< MultiFab > &  S_data,
const MultiFab &  S_prim,
MultiFab &  avg_xmom,
MultiFab &  avg_ymom,
MultiFab &  avg_zmom,
const MultiFab &  xvel,
const MultiFab &  yvel,
const MultiFab &  ,
const MultiFab &  source,
const MultiFab *  SmnSmn,
const 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,
std::unique_ptr< MultiFab > &  z_phys_nd,
std::unique_ptr< MultiFab > &  z_phys_cc,
std::unique_ptr< MultiFab > &  ax,
std::unique_ptr< MultiFab > &  ay,
std::unique_ptr< MultiFab > &  az,
std::unique_ptr< MultiFab > &  detJ,
MultiFab *  detJ_new,
Gpu::DeviceVector< Real > &  stretched_dz_d,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
amrex::EBFArrayBoxFactory const &  ebfact,
YAFluxRegister *  fr_as_crse,
YAFluxRegister *  fr_as_fine 
)

Function for computing the slow RHS for the evolution equations for the scalars other than density or potential temperature

Parameters
[in]evellevel of resolution
[in]finest_levelfinest level of resolution
[in]nrkwhich RK stage
[in]dtslow time step
[out]S_rhsRHS computed here
[in]S_oldsolution at start of time step
[in]S_newsolution at end of current RK stage
[in]S_datacurrent solution
[in]S_primprimitive variables (i.e. conserved variables divided by density)
[in]avg_xmom
[in]avg_ymom
[in]avg_zmom
[in]xvelx-component of velocity
[in]yvely-component of velocity
[in]zvelz-component of velocity
[in]sourcesource terms for conserved variables
[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]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 at start of time step (= 1 if use_terrain is false)
[in]detJ_newJacobian of the metric transformation at new RK stage time (= 1 if use_terrain is false)
[in]mapfacmap factors
[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
101 {
102  BL_PROFILE_REGION("erf_slow_rhs_post()");
103 
104  const BCRec* bc_ptr_d = domain_bcs_type_d.data();
105  const BCRec* bc_ptr_h = domain_bcs_type_h.data();
106 
107  AdvChoice ac = solverChoice.advChoice;
108  DiffChoice dc = solverChoice.diffChoice;
109  TurbChoice tc = solverChoice.turbChoice[level];
110 
111  const MultiFab* t_mean_mf = nullptr;
112  if (SurfLayer) { t_mean_mf = SurfLayer->get_mac_avg(level,2); }
113 
114  const bool l_use_terrain = (solverChoice.mesh_type != MeshType::ConstantDz);
115  const bool l_moving_terrain = (solverChoice.terrain_type == TerrainType::MovingFittedMesh);
116  const bool l_reflux = ( (solverChoice.coupling_type == CouplingType::TwoWay) && (nrk == 2) && (finest_level > 0) );
117  if (l_moving_terrain) AMREX_ALWAYS_ASSERT(l_use_terrain);
118 
119  const bool l_anelastic = solverChoice.anelastic[level];
120 
121  const bool l_use_KE = ( tc.use_tke );
122  const bool l_need_SmnSmn = ( tc.les_type == LESType::Deardorff ||
123  tc.rans_type == RANSType::kEqn );
124  const bool l_advect_KE = ( tc.use_tke && tc.advect_tke );
125  const bool l_use_diff = ((dc.molec_diff_type != MolecDiffType::None) ||
126  (tc.les_type != LESType::None) ||
127  (tc.rans_type != RANSType::None) ||
128  (tc.pbl_type != PBLType::None) );
129  const bool l_use_turb = ( tc.les_type == LESType::Smagorinsky ||
130  tc.les_type == LESType::Deardorff ||
131  tc.rans_type == RANSType::kEqn ||
132  tc.pbl_type == PBLType::MYJ ||
133  tc.pbl_type == PBLType::MYNN25 ||
134  tc.pbl_type == PBLType::MYNNEDMF ||
135  tc.pbl_type == PBLType::YSU ||
136  tc.pbl_type == PBLType::MRF );
137  const bool l_rotate = (solverChoice.use_rotate_surface_flux);
138  const bool do_upwind = solverChoice.upwind_real_bcs;
139  amrex::ignore_unused(do_upwind);
140 
141  const Box& domain = geom.Domain();
142 
143  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
144  const Real* dx = geom.CellSize();
145 
146  // *************************************************************************
147  // Set gravity as a vector
148  // *************************************************************************
149  const Array<Real,AMREX_SPACEDIM> grav{0.0, 0.0, -solverChoice.gravity};
150  const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
151 
152  // *************************************************************************
153  // Pre-computed quantities
154  // *************************************************************************
155  int nvars = S_data[IntVars::cons].nComp();
156  const BoxArray& ba = S_data[IntVars::cons].boxArray();
157  const DistributionMapping& dm = S_data[IntVars::cons].DistributionMap();
158 
159  std::unique_ptr<MultiFab> dflux_x;
160  std::unique_ptr<MultiFab> dflux_y;
161  std::unique_ptr<MultiFab> dflux_z;
162 
163  if (l_use_diff) {
164  dflux_x = std::make_unique<MultiFab>(convert(ba,IntVect(1,0,0)), dm, 1, 0);
165  dflux_y = std::make_unique<MultiFab>(convert(ba,IntVect(0,1,0)), dm, 1, 0);
166  dflux_z = std::make_unique<MultiFab>(convert(ba,IntVect(0,0,1)), dm, 1, 0);
167  } else {
168  dflux_x = nullptr;
169  dflux_y = nullptr;
170  dflux_z = nullptr;
171  }
172 
173  // Valid vars
174  Vector<int> is_valid_slow_var; is_valid_slow_var.resize(RhoQ1_comp+1,0);
175  if (l_use_KE) {is_valid_slow_var[ RhoKE_comp] = 1;}
176  is_valid_slow_var[RhoScalar_comp] = 1;
177  if (solverChoice.moisture_type != MoistureType::None) {
178  is_valid_slow_var[RhoQ1_comp] = 1;
179  }
180 
181  // *************************************************************************
182  // Calculate cell-centered eddy viscosity & diffusivities
183  //
184  // Notes -- we fill all the data in ghost cells before calling this so
185  // that we can fill the eddy viscosity in the ghost regions and
186  // not have to call a boundary filler on this data itself
187  //
188  // LES - updates both horizontal and vertical eddy viscosityS_tmp components
189  // PBL - only updates vertical eddy viscosity components so horizontal
190  // components come from the LES model or are left as zero.
191  // *************************************************************************
192 
193  // *************************************************************************
194  // Define updates and fluxes in the current RK stage
195  // *************************************************************************
196 #ifdef _OPENMP
197 #pragma omp parallel if (Gpu::notInLaunchRegion())
198 #endif
199  {
200  std::array<FArrayBox,AMREX_SPACEDIM> flux;
201 
202  int start_comp;
203  int num_comp;
204 
205  // Cell-centered masks for EB (used for flux interpolation)
206  iMultiFab physbnd_mask;
207  bool already_on_centroids = false;
208  if (solverChoice.terrain_type == TerrainType::EB) {
209  physbnd_mask.define(S_data[IntVars::cons].boxArray(), S_data[IntVars::cons].DistributionMap(), 1, 1);
210  physbnd_mask.BuildMask(geom.Domain(), geom.periodicity(), 1, 1, 0, 1);
211  }
212 
213  for (MFIter mfi(S_data[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
214 
215  Box tbx = mfi.tilebox();
216 
217  // *************************************************************************
218  // Define flux arrays for use in advection
219  // *************************************************************************
220  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
221  if (solverChoice.terrain_type != TerrainType::EB) {
222  flux[dir].resize(surroundingNodes(tbx,dir),nvars,The_Async_Arena());
223  } else {
224  flux[dir].resize(surroundingNodes(tbx,dir).grow(1),nvars,The_Async_Arena());
225  }
226  flux[dir].setVal<RunOn::Device>(0.);
227  }
228  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
229  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
230 
231  // *************************************************************************
232  // Define Array4's
233  // *************************************************************************
234  const Array4<const Real> & old_cons = S_old[IntVars::cons].array(mfi);
235  const Array4< Real> & cell_rhs = S_rhs[IntVars::cons].array(mfi);
236 
237  const Array4< Real> & new_cons = S_new[IntVars::cons].array(mfi);
238  const Array4< Real> & new_xmom = S_new[IntVars::xmom].array(mfi);
239  const Array4< Real> & new_ymom = S_new[IntVars::ymom].array(mfi);
240  const Array4< Real> & new_zmom = S_new[IntVars::zmom].array(mfi);
241 
242  const Array4< Real> & cur_cons = S_data[IntVars::cons].array(mfi);
243  const Array4<const Real> & cur_prim = S_prim.array(mfi);
244  const Array4< Real> & cur_xmom = S_data[IntVars::xmom].array(mfi);
245  const Array4< Real> & cur_ymom = S_data[IntVars::ymom].array(mfi);
246  const Array4< Real> & cur_zmom = S_data[IntVars::zmom].array(mfi);
247 
248  Array4<Real> avg_xmom_arr = avg_xmom.array(mfi);
249  Array4<Real> avg_ymom_arr = avg_ymom.array(mfi);
250  Array4<Real> avg_zmom_arr = avg_zmom.array(mfi);
251 
252  const Array4<const Real> & u = xvel.array(mfi);
253  const Array4<const Real> & v = yvel.array(mfi);
254 
255  const Array4<Real const>& mu_turb = l_use_turb ? eddyDiffs->const_array(mfi) : Array4<const Real>{};
256 
257  const Array4<const Real>& z_nd = z_phys_nd->const_array(mfi);
258  const Array4<const Real>& z_cc = z_phys_cc->const_array(mfi);
259  const Array4<const Real>& detJ_new_arr = l_moving_terrain ? detJ_new->const_array(mfi) : Array4<const Real>{};
260 
261  // Map factors
262  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
263  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
264  const Array4<const Real>& mf_vx = mapfac[MapFacType::v_x]->const_array(mfi);
265  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
266  const Array4<const Real>& mf_uy = mapfac[MapFacType::u_y]->const_array(mfi);
267  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
268 
269  // SmnSmn for KE src with Deardorff or k-eqn RANS
270  const Array4<const Real>& SmnSmn_a = l_need_SmnSmn ? SmnSmn->const_array(mfi) : Array4<const Real>{};
271 
272  // **************************************************************************
273  // Here we fill the "current" data with "new" data because that is the result of the previous RK stage
274  // **************************************************************************
275  int nsv = S_old[IntVars::cons].nComp() - 2;
276  const GpuArray<int, IntVars::NumTypes> scomp_slow = { 2,0,0,0};
277  const GpuArray<int, IntVars::NumTypes> ncomp_slow = {nsv,0,0,0};
278 
279  // **************************************************************************
280  // Note that here we do copy only the "slow" variables, not (rho) or (rho theta)
281  // **************************************************************************
282  ParallelFor(tbx, ncomp_slow[IntVars::cons],
283  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) {
284  const int n = scomp_slow[IntVars::cons] + nn;
285  cur_cons(i,j,k,n) = new_cons(i,j,k,n);
286  });
287 
288  // We have projected the velocities stored in S_data but we will use
289  // the velocities stored in {avg_xmom,avg_ymom,avg_zmom} to update the scalars,
290  // so we need to copy from S_data (projected) into these
291  if (l_anelastic) {
292  Box tbx_inc = mfi.nodaltilebox(0);
293  Box tby_inc = mfi.nodaltilebox(1);
294  Box tbz_inc = mfi.nodaltilebox(2);
295 
296  ParallelFor(tbx_inc, tby_inc, tbz_inc,
297  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
298  avg_xmom_arr(i,j,k) = cur_xmom(i,j,k);
299  },
300  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
301  avg_ymom_arr(i,j,k) = cur_ymom(i,j,k);
302  },
303  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
304  avg_zmom_arr(i,j,k) = cur_zmom(i,j,k);
305  });
306  }
307 
308  // **************************************************************************
309  // Define updates in the RHS of continuity, temperature, and scalar equations
310  // **************************************************************************
311  bool l_eb_terrain_cc = false; // EB terrain on cell-centered grid
312  Array4<const int> mask_arr{};
313  Array4<const EBCellFlag> cfg_arr{};
314  Array4<const Real> ax_arr{};
315  Array4<const Real> ay_arr{};
316  Array4<const Real> az_arr{};
317  Array4<const Real> fcx_arr{};
318  Array4<const Real> fcy_arr{};
319  Array4<const Real> fcz_arr{};
320  Array4<const Real> detJ_arr{};
321  if (solverChoice.terrain_type == TerrainType::EB) {
322  EBCellFlagFab const& cfg = ebfact.getMultiEBCellFlagFab()[mfi];
323  cfg_arr = cfg.const_array();
324  if (cfg.getType(tbx) == FabType::singlevalued) {
325  l_eb_terrain_cc = true;
326  ax_arr = ebfact.getAreaFrac()[0]->const_array(mfi);
327  ay_arr = ebfact.getAreaFrac()[1]->const_array(mfi);
328  az_arr = ebfact.getAreaFrac()[2]->const_array(mfi);
329  fcx_arr = ebfact.getFaceCent()[0]->const_array(mfi);
330  fcy_arr = ebfact.getFaceCent()[1]->const_array(mfi);
331  fcz_arr = ebfact.getFaceCent()[2]->const_array(mfi);
332  detJ_arr = ebfact.getVolFrac().const_array(mfi);
333  // if (!already_on_centroids) {mask_arr = physbnd_mask.const_array(mfi);}
334  mask_arr = physbnd_mask.const_array(mfi);
335  }
336  }
337  if (!l_eb_terrain_cc) {
338  ax_arr = ax->const_array(mfi);
339  ay_arr = ay->const_array(mfi);
340  az_arr = az->const_array(mfi);
341  detJ_arr = detJ->const_array(mfi);
342  }
343 
344  AdvType horiz_adv_type, vert_adv_type;
345  Real horiz_upw_frac, vert_upw_frac;
346 
347  Array4<Real> diffflux_x, diffflux_y, diffflux_z;
348  Array4<Real> hfx_x, hfx_y, hfx_z, diss;
349  Array4<Real> q1fx_x, q1fx_y, q1fx_z, q2fx_z;
350  const bool use_SurfLayer = (SurfLayer != nullptr);
351 
352  if (l_use_diff) {
353  diffflux_x = dflux_x->array(mfi);
354  diffflux_y = dflux_y->array(mfi);
355  diffflux_z = dflux_z->array(mfi);
356 
357  hfx_x = Hfx1->array(mfi);
358  hfx_y = Hfx2->array(mfi);
359  hfx_z = Hfx3->array(mfi);
360  diss = Diss->array(mfi);
361 
362  if (Q1fx1) q1fx_x = Q1fx1->array(mfi);
363  if (Q1fx2) q1fx_y = Q1fx2->array(mfi);
364  if (Q1fx3) q1fx_z = Q1fx3->array(mfi);
365  if (Q2fx3) q2fx_z = Q2fx3->array(mfi);
366  }
367 
368  //
369  // Note that we either advect and diffuse all or none of the moisture variables
370  //
371  for (int ivar(RhoKE_comp); ivar<= RhoQ1_comp; ++ivar)
372  {
373  if (is_valid_slow_var[ivar])
374  {
375  start_comp = ivar;
376  num_comp = 1;
377 
378  if (ivar == RhoQ1_comp) {
379  horiz_adv_type = ac.moistscal_horiz_adv_type;
380  vert_adv_type = ac.moistscal_vert_adv_type;
381  horiz_upw_frac = ac.moistscal_horiz_upw_frac;
382  vert_upw_frac = ac.moistscal_vert_upw_frac;
383 
384  if (ac.use_efficient_advection){
385  horiz_adv_type = EfficientAdvType(nrk,ac.moistscal_horiz_adv_type);
386  vert_adv_type = EfficientAdvType(nrk,ac.moistscal_vert_adv_type);
387  }
388 
389  num_comp = n_qstate;
390 
391  } else {
392  horiz_adv_type = ac.dryscal_horiz_adv_type;
393  vert_adv_type = ac.dryscal_vert_adv_type;
394  horiz_upw_frac = ac.dryscal_horiz_upw_frac;
395  vert_upw_frac = ac.dryscal_vert_upw_frac;
396 
397  if (ac.use_efficient_advection){
398  horiz_adv_type = EfficientAdvType(nrk,ac.dryscal_horiz_adv_type);
399  vert_adv_type = EfficientAdvType(nrk,ac.dryscal_vert_adv_type);
400  }
401 
402  if (ivar == RhoScalar_comp) {
403  num_comp = NSCALARS;
404  }
405  }
406 
407  if (( ivar != RhoKE_comp ) ||
408  ((ivar == RhoKE_comp) && l_advect_KE))
409  {
410  if (!l_eb_terrain_cc){
411  AdvectionSrcForScalars(tbx, start_comp, num_comp,
412  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
413  cur_prim, cell_rhs,
414  detJ_arr, dxInv, mf_mx, mf_my,
415  horiz_adv_type, vert_adv_type,
416  horiz_upw_frac, vert_upw_frac,
417  flx_arr, domain, bc_ptr_h);
418  } else {
419  EBAdvectionSrcForScalars(tbx, start_comp, num_comp,
420  avg_xmom_arr, avg_ymom_arr, avg_zmom_arr,
421  cur_prim, cell_rhs,
422  mask_arr, cfg_arr, ax_arr, ay_arr, az_arr,
423  fcx_arr, fcy_arr, fcz_arr,
424  detJ_arr, dxInv, mf_mx, mf_my,
425  horiz_adv_type, vert_adv_type,
426  horiz_upw_frac, vert_upw_frac,
427  flx_arr, domain, bc_ptr_h,
428  already_on_centroids);
429  }
430  }
431 
432  if (l_use_diff)
433  {
434  // Here we hardwire this to 0 because we only use vert_implicit_fac for (rho_theta)
435  const Real l_vert_implicit_fac = 0.0;;
436 
437  const Array4<const Real> tm_arr = t_mean_mf ? t_mean_mf->const_array(mfi) : Array4<const Real>{};
438 
439  if (solverChoice.mesh_type == MeshType::StretchedDz && solverChoice.terrain_type != TerrainType::EB) {
440  DiffusionSrcForState_S(tbx, domain, start_comp, num_comp, u, v,
441  new_cons, cur_prim, cell_rhs,
442  diffflux_x, diffflux_y, diffflux_z,
443  stretched_dz_d, dxInv, SmnSmn_a,
444  mf_mx, mf_ux, mf_vx,
445  mf_my, mf_uy, mf_vy,
446  hfx_z, q1fx_z, q2fx_z, diss,
447  mu_turb, solverChoice, level,
448  tm_arr, grav_gpu, bc_ptr_d, use_SurfLayer, l_vert_implicit_fac);
449  } else if (l_use_terrain) {
450  DiffusionSrcForState_T(tbx, domain, start_comp, num_comp, l_rotate, u, v,
451  new_cons, cur_prim, cell_rhs,
452  diffflux_x, diffflux_y, diffflux_z,
453  z_nd, z_cc, ax_arr, ay_arr, az_arr,
454  detJ_arr, dxInv, SmnSmn_a,
455  mf_mx, mf_ux, mf_vx,
456  mf_my, mf_uy, mf_vy,
457  hfx_x, hfx_y, hfx_z, q1fx_x, q1fx_y, q1fx_z,q2fx_z, diss,
458  mu_turb, solverChoice, level,
459  tm_arr, grav_gpu, bc_ptr_d, use_SurfLayer, l_vert_implicit_fac);
460  } else {
461  DiffusionSrcForState_N(tbx, domain, start_comp, num_comp, u, v,
462  new_cons, cur_prim, cell_rhs,
463  diffflux_x, diffflux_y, diffflux_z, dxInv, SmnSmn_a,
464  mf_mx, mf_ux, mf_vx,
465  mf_my, mf_uy, mf_vy,
466  hfx_z, q1fx_z, q2fx_z, diss,
467  mu_turb, solverChoice, level,
468  tm_arr, grav_gpu, bc_ptr_d, use_SurfLayer, l_vert_implicit_fac);
469  }
470  } // use_diff
471  } // valid slow var
472  } // loop ivar
473 
474 #if defined(ERF_USE_NETCDF)
475  if (moist_set_rhs_bool)
476  {
477  const Array4<const Real> & old_cons_const = S_old[IntVars::cons].const_array(mfi);
478  const Array4<const Real> & new_cons_const = S_new[IntVars::cons].const_array(mfi);
479  moist_set_rhs(geom, tbx, old_cons_const, new_cons_const, cell_rhs, bdy_time_interval,
480  new_stage_time, dt, stop_time_elapsed, width, set_width, do_upwind, domain,
481  bdy_data_xlo, bdy_data_xhi, bdy_data_ylo, bdy_data_yhi);
482  }
483 #endif
484 
485 #ifdef ERF_USE_SHOC
486  if (solverChoice.use_shoc) {
487  shoc_lev->add_slow_tend(mfi,tbx,cell_rhs);
488  }
489 #endif
490 
491  // This updates just the "slow" conserved variables
492  {
493  BL_PROFILE("rhs_post_8");
494 
496 
497  auto const& src_arr = source.const_array(mfi);
498 
499  for (int ivar(RhoKE_comp); ivar<= RhoQ1_comp; ++ivar)
500  {
501  if (is_valid_slow_var[ivar])
502  {
503  start_comp = ivar;
504  num_comp = 1;
505  if (ivar == RhoQ1_comp) {
506  num_comp = nvars - RhoQ1_comp;
507  } else if (ivar == RhoScalar_comp) {
508  num_comp = NSCALARS;
509  }
510 
511  if (l_moving_terrain)
512  {
513  ParallelFor(tbx, num_comp,
514  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) noexcept {
515  const int n = start_comp + nn;
516  cell_rhs(i,j,k,n) += src_arr(i,j,k,n);
517  Real temp_val = detJ_arr(i,j,k) * old_cons(i,j,k,n) + dt * detJ_arr(i,j,k) * cell_rhs(i,j,k,n);
518  cur_cons(i,j,k,n) = temp_val / detJ_new_arr(i,j,k);
519  if (ivar == RhoKE_comp) {
520  cur_cons(i,j,k,n) = amrex::max(cur_cons(i,j,k,n), eps);
521  }
522  });
523 
524  } else if (l_anelastic && (nrk == 1)) { // not moving and ( (anelastic) and second RK stage) )
525 
526  ParallelFor(tbx, num_comp,
527  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) noexcept {
528  const int n = start_comp + nn;
529  cell_rhs(i,j,k,n) += src_arr(i,j,k,n);
530 
531  // Re-construct the cell_rhs used in the first RK stage
532  Real dt_times_old_cell_rhs = cur_cons(i,j,k,n) - old_cons(i,j,k,n);
533 
534  // Add the time-averaged RHS to the old state
535  cur_cons(i,j,k,n) = old_cons(i,j,k,n) + 0.5 * (dt_times_old_cell_rhs + dt * cell_rhs(i,j,k,n));
536 
537  if (ivar == RhoKE_comp) {
538  cur_cons(i,j,k,n) = amrex::max(cur_cons(i,j,k,n), eps);
539  } else if (ivar >= RhoQ1_comp) {
540  cur_cons(i,j,k,n) = amrex::max(cur_cons(i,j,k,n), 0.0);
541  }
542  });
543 
544  } else { // not moving and ( (not anelastic) or (first RK stage) )
545 
546  ParallelFor(tbx, num_comp,
547  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) noexcept {
548  const int n = start_comp + nn;
549  cell_rhs(i,j,k,n) += src_arr(i,j,k,n);
550  cur_cons(i,j,k,n) = old_cons(i,j,k,n) + dt * cell_rhs(i,j,k,n);
551  if (ivar == RhoKE_comp) {
552  cur_cons(i,j,k,n) = amrex::max(cur_cons(i,j,k,n), eps);
553  } else if (ivar >= RhoQ1_comp) {
554  cur_cons(i,j,k,n) = amrex::max(cur_cons(i,j,k,n), 0.0);
555  }
556  });
557 
558  } // moving, anelastic or neither?
559 
560  } // is_valid
561  } // ivar
562  } // profile
563 
564  {
565  BL_PROFILE("rhs_post_9");
566  // This updates all the conserved variables (not just the "slow" ones)
567  int num_comp_all = S_data[IntVars::cons].nComp();
568  ParallelFor(tbx, num_comp_all,
569  [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept {
570  new_cons(i,j,k,n) = cur_cons(i,j,k,n);
571  });
572  } // end profile
573 
574  Box xtbx = mfi.nodaltilebox(0);
575  Box ytbx = mfi.nodaltilebox(1);
576  Box ztbx = mfi.nodaltilebox(2);
577 
578  {
579  BL_PROFILE("rhs_post_10()");
580  ParallelFor(xtbx, ytbx, ztbx,
581  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
582  new_xmom(i,j,k) = cur_xmom(i,j,k);
583  },
584  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
585  new_ymom(i,j,k) = cur_ymom(i,j,k);
586  },
587  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
588  new_zmom(i,j,k) = cur_zmom(i,j,k);
589  });
590  } // end profile
591 
592  {
593  BL_PROFILE("rhs_post_10");
594  // We only add to the flux registers in the final RK step
595  if (l_reflux) {
596  int strt_comp_reflux = RhoTheta_comp + 1;
597  int num_comp_reflux = nvars - strt_comp_reflux;
598  if (level < finest_level) {
599  fr_as_crse->CrseAdd(mfi,
600  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
601  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
602  }
603  if (level > 0) {
604  fr_as_fine->FineAdd(mfi,
605  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
606  dx, dt, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
607  }
608 
609  // This is necessary here so we don't go on to the next FArrayBox without
610  // having finished copying the fluxes into the FluxRegisters (since the fluxes
611  // are stored in temporary FArrayBox's)
612  Gpu::streamSynchronize();
613 
614  } // two-way coupling
615  } // end profile
616  } // mfi
617  } // OMP
618 }
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)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE AdvType EfficientAdvType(int nrk, AdvType adv_type)
Definition: ERF_Advection.H:281
@ nvars
Definition: ERF_DataStruct.H:97
@ v_x
Definition: ERF_DataStruct.H:24
@ u_y
Definition: ERF_DataStruct.H:25
@ v_y
Definition: ERF_DataStruct.H:25
@ m_y
Definition: ERF_DataStruct.H:25
@ u_x
Definition: ERF_DataStruct.H:24
@ m_x
Definition: ERF_DataStruct.H:24
void 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)
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
AdvType
Definition: ERF_IndexDefines.H:221
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
@ xvel
Definition: ERF_IndexDefines.H:141
@ yvel
Definition: ERF_IndexDefines.H:142
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
real(c_double), private ac
Definition: ERF_module_mp_morr_two_moment.F90:181
Definition: ERF_AdvStruct.H:19
Definition: ERF_DiffStruct.H:19
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
static MeshType mesh_type
Definition: ERF_DataStruct.H:1048
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1057
amrex::Real gravity
Definition: ERF_DataStruct.H:1122
bool use_shoc
Definition: ERF_DataStruct.H:1155
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1060
amrex::Vector< int > anelastic
Definition: ERF_DataStruct.H:1066
static bool upwind_real_bcs
Definition: ERF_DataStruct.H:1045
AdvChoice advChoice
Definition: ERF_DataStruct.H:1056
MoistureType moisture_type
Definition: ERF_DataStruct.H:1171
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1036
bool use_rotate_surface_flux
Definition: ERF_DataStruct.H:1152
CouplingType coupling_type
Definition: ERF_DataStruct.H:1170
Definition: ERF_TurbStruct.H:42
PBLType pbl_type
Definition: ERF_TurbStruct.H:418
RANSType rans_type
Definition: ERF_TurbStruct.H:413
bool advect_tke
Definition: ERF_TurbStruct.H:463
bool use_tke
Definition: ERF_TurbStruct.H:430
LESType les_type
Definition: ERF_TurbStruct.H:371
Here is the call graph for this function: