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