ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TI_slow_rhs_pre.H
Go to the documentation of this file.
1 #include "ERF_SrcHeaders.H"
2 
3 /**
4  * Wrapper for calling the routine that creates the slow RHS
5  */
6  auto slow_rhs_fun_pre = [&](Vector<MultiFab>& S_rhs,
7  Vector<MultiFab>& S_old,
8  Vector<MultiFab>& S_data,
9  const Real old_step_time,
10  const Real old_stage_time,
11  const Real new_stage_time,
12  const int nrk)
13  {
14  BL_PROFILE("slow_rhs_fun_pre");
15  //
16  // Define primitive variables for all later RK stages
17  // (We have already done this for the first RK step)
18  // Note that it is essential this happen before the call to make_mom_sources
19  // because some of the buoyancy routines use the primitive variables
20  //
21  if (nrk > 0) {
22  int ng_cons = S_data[IntVars::cons].nGrow();
23  cons_to_prim(S_data[IntVars::cons], S_prim, ng_cons);
25  }
26 
27  if (verbose) Print() << std::setprecision(timeprecision)
28  << "Making slow rhs at time " << old_stage_time
29  << " for fast variables advancing from " << old_step_time
30  << " to " << new_stage_time << std::endl;
31 
32  Real slow_dt = new_stage_time - old_step_time;
33 
34  const GpuArray<Real, AMREX_SPACEDIM> dxInv = fine_geom.InvCellSizeArray();
35 
36  // *************************************************************************
37  // Set up flux registers if using two_way coupling
38  // *************************************************************************
39  YAFluxRegister* fr_as_crse = nullptr;
40  YAFluxRegister* fr_as_fine = nullptr;
41  if (solverChoice.coupling_type == CouplingType::TwoWay && finest_level > 0) {
42  if (level < finest_level) {
43  fr_as_crse = getAdvFluxReg(level+1);
44  fr_as_crse->reset();
45  }
46  if (level > 0) {
47  fr_as_fine = getAdvFluxReg(level);
48  }
49  }
50 
51  // *************************************************************************
52  // Get multifab pointers
53  // *************************************************************************
54 
55  // Canopy data for mom sources
56  MultiFab* forest_drag = (solverChoice.do_forest_drag) ?
57  m_forest_drag[level]->get_drag_field() : nullptr;
58 
59  // Immersed Forcing
60  MultiFab* terrain_blank = (solverChoice.terrain_type == TerrainType::ImmersedForcing ||
61  solverChoice.buildings_type == BuildingsType::ImmersedForcing) ?
62  terrain_blanking[level].get() : nullptr;
63 
64  // Update the total moisture variable *before* computing sources since this is used in
65  // the buoyancy calculation
66  if (solverChoice.moisture_type != MoistureType::None) {
67  int n_qstate_into_total = micro->Get_Qstate_Moist_Size() - micro->Get_Qstate_Moist_NumConc_Size();
68  make_qt(S_data[IntVars::cons], qt, n_qstate_into_total);
69  }
70 
71  MultiFab p0_to_use, base_to_use;
72  MultiFab *zpn_to_use, *zpc_to_use, *ax_to_use, *ay_to_use, *az_to_use, *dJ_to_use;
73 
74  // Moving terrain
75  std::unique_ptr<MultiFab> z_t_pert;
76  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
77  {
78  z_t_pert = std::make_unique<MultiFab>(S_data[IntVars::zmom].boxArray(), S_data[IntVars::zmom].DistributionMap(), 1, 1);
79  update_terrain_stage(level, old_step_time, old_stage_time, new_stage_time, slow_dt);
80 
81  p0_to_use = MultiFab(base_state_new[level], make_alias, BaseState::p0_comp, 1);
82  base_to_use = MultiFab(base_state_new[level], make_alias, 0, BaseState::num_comps);
83  zpn_to_use = z_phys_nd_src[level].get();
84  zpc_to_use = z_phys_cc_src[level].get();
85  ax_to_use = ax_src[level].get();
86  ay_to_use = ay_src[level].get();
87  az_to_use = az_src[level].get();
88  dJ_to_use = detJ_cc_src[level].get();
89 
90  } else {
91  p0_to_use = MultiFab(base_state[level], make_alias, BaseState::p0_comp, 1);
92  base_to_use = MultiFab(base_state[level], make_alias, 0, BaseState::num_comps);
93  zpn_to_use = z_phys_nd[level].get();
94  zpc_to_use = z_phys_cc[level].get();
95  ax_to_use = ax[level].get();
96  ay_to_use = ay[level].get();
97  az_to_use = az[level].get();
98  dJ_to_use = detJ_cc[level].get();
99  }
100 
101  // *****************************************************************************
102  // Construct the source terms for the cell-centered (conserved) variables
103  // *****************************************************************************
104  make_sources(level, nrk, slow_dt, old_stage_time,
105  S_data, S_prim, cc_src, base_state[level], zpc_to_use,
106  xvel_new, yvel_new,
107  qheating_rates[level].get(),
108  terrain_blank, fine_geom, solverChoice,
109  mapfac[level],
110  rhotheta_src[level].get(), rhoqt_src[level].get(),
111  dptr_wbar_sub, d_rayleigh_ptrs_at_lev, d_sinesq_at_lev,
112  (solverChoice.hindcast_surface_bcs? &surface_state_interp[level] : nullptr),
113  input_sounding_data, turbPert, true);
114 
115  // *****************************************************************************
116  // Define the pressure gradient
117  // *****************************************************************************
118  make_gradp_pert(level, solverChoice, fine_geom, S_data,
119  p0_to_use, *zpn_to_use, *zpc_to_use,
120  mapfac[level],
121  get_eb(level), gradp[level]);
122 
123  // *****************************************************************************
124  // Define the buoyancy forcing term in the z-direction
125  // *****************************************************************************
126  int num_q = micro->Get_Qstate_Moist_Size() - micro->Get_Qstate_Moist_NumConc_Size();
127  make_buoyancy(level, S_data, S_prim, qt, buoyancy, fine_geom, solverChoice, base_to_use,
128  num_q, get_eb(level), solverChoice.anelastic[level]);
129 
130  // *****************************************************************************
131  // Make remaining (not gradp or buoyancy) momentum sources
132  // *****************************************************************************
133  make_mom_sources(old_stage_time, slow_dt,
134  S_data, zpn_to_use, zpc_to_use, stretched_dz_h[level],
135  xvel_new, yvel_new, zvel_new,
136  xmom_src, ymom_src, zmom_src,
137  base_to_use, forest_drag, terrain_blank,
138  cosPhi_m[level].get(), sinPhi_m[level].get(), fine_geom, solverChoice,
139  mapfac[level],
140  (solverChoice.have_geo_wind_profile) ? d_u_geos[level].data(): nullptr,
141  (solverChoice.have_geo_wind_profile) ? d_v_geos[level].data(): nullptr,
142  dptr_wbar_sub, d_rayleigh_ptrs_at_lev, d_sinesq_at_lev, d_sinesq_stag_at_lev,
143  d_sponge_ptrs_at_lev,
144  (solverChoice.hindcast_lateral_forcing? &forecast_state_interp[level] : nullptr),
145  (solverChoice.hindcast_surface_bcs? &surface_state_interp[level] : nullptr),
146  input_sounding_data, get_eb(level), true);
147 
148  // *****************************************************************************
149  // Add body sources if doing flow around a body
150  // *****************************************************************************
151  add_thin_body_sources(xmom_src, ymom_src, zmom_src,
152  xflux_imask[level], yflux_imask[level], zflux_imask[level],
153  thin_xforce[level], thin_yforce[level], thin_zforce[level]);
154 
155  // *****************************************************************************
156  // Define RHS for rho, rho_theta and momenta
157  // *****************************************************************************
158  erf_slow_rhs_pre(level, finest_level, nrk, slow_dt, S_rhs, S_old, S_data,
159  S_prim, qt, avg_xmom[level], avg_ymom[level], avg_zmom[level],
160  xvel_new, yvel_new, zvel_new,
161  z_t_rk[level], cc_src, xmom_src, ymom_src, zmom_src, buoyancy,
162  (level > 0) ? &zmom_crse_rhs[level] : nullptr,
163  Tau[level], Tau_corr[level], Tau_EB[level],
164  SmnSmn, eddyDiffs, Hfx1, Hfx2, Hfx3, Q1fx1, Q1fx2, Q1fx3, Q2fx3, Diss, Hfx3_EB,
165  fine_geom, solverChoice, m_SurfaceLayer, domain_bcs_type_d, domain_bcs_type,
166  *zpn_to_use, *zpc_to_use, *ax_to_use, *ay_to_use, *az_to_use, *dJ_to_use,
167  stretched_dz_d[level], gradp[level],
168  mapfac[level], get_eb(level),
169 #ifdef ERF_USE_SHOC
170  shoc_interface[level],
171 #endif
172  fr_as_crse, fr_as_fine);
173 
174  if ((solverChoice.vert_implicit_fac[nrk] > zero) && solverChoice.implicit_before_substep) {
175  MultiFab scratch(S_data[IntVars::cons].boxArray(),S_data[IntVars::cons].DistributionMap(), 2,
176  S_data[IntVars::cons].nGrowVect());
177  MultiFab::Copy(scratch, S_old[IntVars::cons], 0, 0, 2, S_data[IntVars::cons].nGrowVect()); // scratch := S_old (for rho, rhotheta)
178  MultiFab::Saxpy(scratch, slow_dt, S_rhs[IntVars::cons], 0, 0, 2, 0); // scratch := S_old + slow_dt*Src (for rho, rhotheta)
179  scratch.FillBoundary(geom[level].periodicity());
180 
181  MultiFab scratch_xmom(S_data[IntVars::xmom].boxArray(),
182  S_data[IntVars::xmom].DistributionMap(), 1,
183  S_data[IntVars::xmom].nGrowVect());
184  MultiFab scratch_ymom(S_data[IntVars::ymom].boxArray(),
185  S_data[IntVars::ymom].DistributionMap(), 1,
186  S_data[IntVars::ymom].nGrowVect());
187 #ifdef ERF_IMPLICIT_W
188  MultiFab scratch_zmom(S_data[IntVars::zmom].boxArray(),
189  S_data[IntVars::zmom].DistributionMap(), 1,
190  S_data[IntVars::zmom].nGrowVect());
191 #endif
192  if (solverChoice.implicit_momentum_diffusion) {
193  MultiFab::Copy(scratch_xmom, S_old[IntVars::xmom], 0, 0, 1, S_data[IntVars::xmom].nGrowVect()); // scratch := S_old
194  MultiFab::Saxpy(scratch_xmom, slow_dt, S_rhs[IntVars::xmom], 0, 0, 1, 0); // scratch := S_old + slow_dt*Src
195  scratch_xmom.FillBoundary(geom[level].periodicity());
196 
197  MultiFab::Copy(scratch_ymom, S_old[IntVars::ymom], 0, 0, 1, S_data[IntVars::ymom].nGrowVect()); // scratch := S_old
198  MultiFab::Saxpy(scratch_ymom, slow_dt, S_rhs[IntVars::ymom], 0, 0, 1, 0); // scratch := S_old + slow_dt*Src
199  scratch_ymom.FillBoundary(geom[level].periodicity());
200 #ifdef ERF_IMPLICIT_W
201  MultiFab::Copy(scratch_zmom, S_old[IntVars::zmom], 0, 0, 1, S_data[IntVars::zmom].nGrowVect()); // scratch := S_old
202  MultiFab::Saxpy(scratch_zmom, slow_dt, S_rhs[IntVars::zmom], 0, 0, 1, 0); // scratch := S_old + slow_dt*Src
203  scratch_zmom.FillBoundary(geom[level].periodicity());
204 #endif
205  }
206 
207 #include "ERF_ImplicitPre.H"
208 
209  MultiFab::Saxpy(scratch, -one, S_old[IntVars::cons], 1, 1, 1, 0); // scratch := (S_new - S_old) (for rhotheta only)
210  scratch.mult(one / slow_dt); // scratch := (S_new - S_old) / slow_dt
211  MultiFab::Copy(S_rhs[IntVars::cons], scratch, 1, 1, 1, 0); // slow_rhs := (S_new - S_old) / slow_dt (for rhotheta only)
212 
213  if (solverChoice.implicit_momentum_diffusion) {
214  MultiFab::Saxpy(scratch_xmom, -one, S_old[IntVars::xmom], 0, 0, 1, 0); // scratch := (S_new - S_old)
215  scratch_xmom.mult(one / slow_dt); // scratch := (S_new - S_old) / slow_dt
216  MultiFab::Copy(S_rhs[IntVars::xmom], scratch_xmom, 0, 0, 1, 0); // slow_rhs := (S_new - S_old) / slow_dt
217 
218  MultiFab::Saxpy(scratch_ymom, -one, S_old[IntVars::ymom], 0, 0, 1, 0); // scratch := (S_new - S_old)
219  scratch_ymom.mult(one / slow_dt); // scratch := (S_new - S_old) / slow_dt
220  MultiFab::Copy(S_rhs[IntVars::ymom], scratch_ymom, 0, 0, 1, 0); // slow_rhs := (S_new - S_old) / slow_dt
221 #ifdef ERF_IMPLICIT_W
222  MultiFab::Saxpy(scratch_zmom, -one, S_old[IntVars::zmom], 0, 0, 1, 0); // scratch := (S_new - S_old)
223  scratch_zmom.mult(one / slow_dt); // scratch := (S_new - S_old) / slow_dt
224  MultiFab::Copy(S_rhs[IntVars::zmom], scratch_zmom, 0, 0, 1, 0); // slow_rhs := (S_new - S_old) / slow_dt
225 #endif
226  }
227  }
228 
229 #ifdef ERF_USE_SHOC
230  if (solverChoice.use_shoc) {
231  shoc_interface[level]->add_fast_tend(S_rhs);
232  }
233 #endif
234 
235  // *****************************************************************************
236  // Update for moving terrain
237  // *****************************************************************************
238  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
239  {
240  MultiFab r_hse_new (base_state_new[level], make_alias, BaseState::r0_comp, 1);
241  MultiFab p_hse_new (base_state_new[level], make_alias, BaseState::p0_comp, 1);
242  MultiFab pi_hse_new (base_state_new[level], make_alias, BaseState::pi0_comp, 1);
243  MultiFab th_hse_new (base_state_new[level], make_alias, BaseState::th0_comp, 1);
244 
245  MultiFab* r0_new = &r_hse_new;
246  MultiFab* p0_new = &p_hse_new;
247  MultiFab* pi0_new = &pi_hse_new;
248  MultiFab* th0_new = &th_hse_new;
249 
250  // We define and evolve (rho theta)_0 in order to re-create p_0 in a way that is consistent
251  // with our update of (rho theta) but does NOT maintain dp_0 / dz = -rho_0 g. This is why
252  // we no longer discretize the vertical pressure gradient in perturbational form.
253  MultiFab rt0(p0->boxArray(),p0->DistributionMap(),1,1);
254  MultiFab rt0_new(p0->boxArray(),p0->DistributionMap(),1,1);
255  MultiFab r0_temp(p0->boxArray(),p0->DistributionMap(),1,1);
256 
257  // Remember this does NOT maintain dp_0 / dz = -rho_0 g, so we can no longer
258  // discretize the vertical pressure gradient in perturbational form.
259  AMREX_ALWAYS_ASSERT(solverChoice.advChoice.dycore_horiz_adv_type == AdvType::Centered_2nd);
260  AMREX_ALWAYS_ASSERT(solverChoice.advChoice.dycore_vert_adv_type == AdvType::Centered_2nd);
261 
262  Real dt_base = (new_stage_time - old_step_time);
263 
264  const Real l_rdOcp = solverChoice.rdOcp;
265 
266 #ifdef _OPENMP
267 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
268 #endif
269  for ( MFIter mfi(*p0,TilingIfNotGPU()); mfi.isValid(); ++mfi)
270  {
271  const Array4<Real > rt0_arr = rt0.array(mfi);
272  const Array4<Real > rt0_tmp_arr = rt0_new.array(mfi);
273 
274  const Array4<Real const> r0_arr = r0->const_array(mfi);
275  const Array4<Real > r0_new_arr = r0_new->array(mfi);
276  const Array4<Real > r0_tmp_arr = r0_temp.array(mfi);
277 
278  const Array4<Real const> p0_arr = p0->const_array(mfi);
279  const Array4<Real > p0_new_arr = p0_new->array(mfi);
280  const Array4<Real > pi0_new_arr = pi0_new->array(mfi);
281  const Array4<Real > th0_new_arr = th0_new->array(mfi);
282 
283  const Array4<Real >& z_t_arr = z_t_rk[level]->array(mfi);
284 
285  const Array4<Real const>& dJ_old_arr = detJ_cc[level]->const_array(mfi);
286  const Array4<Real const>& dJ_new_arr = detJ_cc_new[level]->const_array(mfi);
287  const Array4<Real const>& dJ_src_arr = detJ_cc_src[level]->const_array(mfi);
288 
289  Box gbx = mfi.growntilebox({1,1,1});
290  amrex::ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
291  {
292  rt0_arr(i,j,k) = getRhoThetagivenP(p0_arr(i,j,k));
293  rt0_tmp_arr(i,j,k) = getRhoThetagivenP(p0_new_arr(i,j,k));
294  r0_tmp_arr(i,j,k) = r0_new_arr(i,j,k);
295  });
296 
297  Box gbx2 = mfi.growntilebox({1,1,0});
298  amrex::ParallelFor(gbx2, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
299  {
300  Real zflux_r_lo = -z_t_arr(i,j,k ) * myhalf * (r0_tmp_arr(i,j,k) + r0_tmp_arr(i,j,k-1));
301  Real zflux_r_hi = -z_t_arr(i,j,k+1) * myhalf * (r0_tmp_arr(i,j,k) + r0_tmp_arr(i,j,k+1));
302 
303  Real zflux_rt_lo = zflux_r_lo * myhalf * (rt0_tmp_arr(i,j,k)/r0_tmp_arr(i,j,k) + rt0_tmp_arr(i,j,k-1)/r0_tmp_arr(i,j,k-1));
304  Real zflux_rt_hi = zflux_r_hi * myhalf * (rt0_tmp_arr(i,j,k)/r0_tmp_arr(i,j,k) + rt0_tmp_arr(i,j,k+1)/r0_tmp_arr(i,j,k+1));
305 
306  Real invdetJ = one / dJ_src_arr(i,j,k);
307 
308  Real src_r = - invdetJ * ( zflux_r_hi - zflux_r_lo ) * dxInv[2];
309  Real src_rt = - invdetJ * ( zflux_rt_hi - zflux_rt_lo ) * dxInv[2];
310 
311  Real rho0_new = dJ_old_arr(i,j,k) * r0_arr(i,j,k) + dt_base * dJ_src_arr(i,j,k) * src_r;
312  Real rt0_tmp_new = dJ_old_arr(i,j,k) * rt0_arr(i,j,k) + dt_base * dJ_src_arr(i,j,k) * src_rt;
313 
314  r0_new_arr(i,j,k) = rho0_new / dJ_new_arr(i,j,k);
315  rt0_tmp_new /= dJ_new_arr(i,j,k);
316 
317  p0_new_arr(i,j,k) = getPgivenRTh(rt0_tmp_new);
318  pi0_new_arr(i,j,k) = getExnergivenRTh(rt0_tmp_new, l_rdOcp);
319  th0_new_arr(i,j,k) = rt0_tmp_new / r0_new_arr(i,j,k);
320  });
321  } // MFIter
322  r0_new->FillBoundary(fine_geom.periodicity());
323  p0_new->FillBoundary(fine_geom.periodicity());
324  th0_new->FillBoundary(fine_geom.periodicity());
325  }
326 
327 #ifdef ERF_USE_NETCDF
328  // Populate RHS for relaxation zones if using real bcs
329  if (solverChoice.use_real_bcs && (level == 0)) {
330  const Real bdy_factor = solverChoice.bdy_nudge_factor;
331  if (real_width>0) {
332  //
333  // Note that old_stage_time is elapsed time, but (start_time+old_stage_time) is total time
334  // start_bdy_time and final_bdy_time are total time
335  //
336  Real total_time = start_time + old_stage_time;
337  realbdy_compute_interior_ghost_rhs(total_time, slow_dt,
338  start_bdy_time, final_bdy_time, bdy_time_interval,
339  bdy_factor, real_width, fine_geom,
340  S_rhs, S_data,
341  bdy_data_xlo, bdy_data_xhi,
342  bdy_data_ylo, bdy_data_yhi,
343  m_r2d);
344  }
345  }
346 #endif
347  }; // end slow_rhs_fun_pre
void add_thin_body_sources(MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, std::unique_ptr< iMultiFab > &xflux_imask_lev, std::unique_ptr< iMultiFab > &yflux_imask_lev, std::unique_ptr< iMultiFab > &zflux_imask_lev, std::unique_ptr< MultiFab > &thin_xforce_lev, std::unique_ptr< MultiFab > &thin_yforce_lev, std::unique_ptr< MultiFab > &thin_zforce_lev)
Definition: ERF_AddThinBodySources.cpp:27
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhoThetagivenP(const amrex::Real p, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:172
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:156
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
@ Centered_2nd
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
pp get("wavelength", wavelength)
void realbdy_compute_interior_ghost_rhs(const Real &time, const Real &delta_t, const Real &start_bdy_time, const Real &final_bdy_time, const Real &bdy_time_interval, const Real &nudge_factor, int width, const Geometry &geom, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_cur_data, Vector< Vector< FArrayBox >> &bdy_data_xlo, Vector< Vector< FArrayBox >> &bdy_data_xhi, Vector< Vector< FArrayBox >> &bdy_data_ylo, Vector< Vector< FArrayBox >> &bdy_data_yhi, std::unique_ptr< ReadBndryPlanes > &m_r2d)
Definition: ERF_InteriorGhostCells.cpp:161
void make_buoyancy(int lev, const Vector< MultiFab > &S_data, const MultiFab &S_prim, const MultiFab &qt, MultiFab &buoyancy, const Geometry geom, const SolverChoice &solverChoice, const MultiFab &base_state, const int n_qstate, const eb_ &ebfact, const int anelastic)
Definition: ERF_MakeBuoyancy.cpp:32
void make_gradp_pert(int level, const SolverChoice &solverChoice, const Geometry &geom, Vector< MultiFab > &S_data, const MultiFab &p0, const MultiFab &z_phys_nd, const MultiFab &z_phys_cc, Vector< std::unique_ptr< MultiFab >> &mapfac, const eb_ &ebfact, Vector< MultiFab > &gradp)
Definition: ERF_MakeGradP.cpp:28
void make_mom_sources(Real time, Real, const Vector< MultiFab > &S_data, const MultiFab *z_phys_nd, const MultiFab *z_phys_cc, Vector< Real > &stretched_dz_h, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &wvel, MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, const MultiFab &base_state, MultiFab *forest_drag, MultiFab *terrain_blank, MultiFab *cosPhi_mf, MultiFab *sinPhi_mf, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &, const Real *dptr_u_geos, const Real *dptr_v_geos, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, const amrex::Real *d_sinesq_at_lev, const amrex::Real *d_sinesq_stag_at_lev, const Vector< Real * > d_sponge_ptrs_at_lev, const Vector< MultiFab > *forecast_state_at_lev, const MultiFab *surface_state_at_lev, InputSoundingData &input_sounding_data, const eb_ &ebfact, bool is_slow_step)
Definition: ERF_MakeMomSources.cpp:37
void make_sources(int level, int, Real dt, Real time, const Vector< MultiFab > &S_data, const MultiFab &S_prim, MultiFab &source, const MultiFab &base_state, const MultiFab *z_phys_cc, const MultiFab &xvel, const MultiFab &yvel, const MultiFab *qheating_rates, MultiFab *terrain_blank, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &mapfac, const MultiFab *rhotheta_src, const MultiFab *rhoqt_src, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, const Real *d_sinesq_at_lev, const MultiFab *surface_state_at_lev, InputSoundingData &input_sounding_data, TurbulentPerturbation &turbPert, bool is_slow_step)
Definition: ERF_MakeSources.cpp:34
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
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, Vector< std::unique_ptr< MultiFab >> &Tau_corr_lev, Vector< Vector< std::unique_ptr< MultiFab >>> &Tau_EB, MultiFab *SmnSmn, MultiFab *eddyDiffs, MultiFab *Hfx1, MultiFab *Hfx2, MultiFab *Hfx3, MultiFab *Q1fx1, MultiFab *Q1fx2, MultiFab *Q1fx3, MultiFab *Q2fx3, MultiFab *Diss, MultiFab *Hfx3_EB, 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)
Definition: ERF_SlowRhsPre.cpp:65
auto slow_rhs_fun_pre
Definition: ERF_TI_slow_rhs_pre.H:6
auto make_pi_stage
Definition: ERF_TI_utils.H:4
auto update_terrain_stage
Definition: ERF_TI_utils.H:125
void cons_to_prim(const MultiFab &cons_state, MultiFab &S_prim, int ng)
Definition: ERF_Utils.cpp:6
void make_qt(const MultiFab &cons_state, MultiFab &qt, int n_qstate_into_total)
Definition: ERF_Utils.cpp:37
@ num_comps
Definition: ERF_IndexDefines.H:78
@ pi0_comp
Definition: ERF_IndexDefines.H:75
@ p0_comp
Definition: ERF_IndexDefines.H:74
@ th0_comp
Definition: ERF_IndexDefines.H:76
@ r0_comp
Definition: ERF_IndexDefines.H:73
@ ymom
Definition: ERF_IndexDefines.H:194
@ cons
Definition: ERF_IndexDefines.H:192
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193
@ qt
Definition: ERF_Kessler.H:28
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40
real(kind=kind_phys), parameter, private r0
Definition: ERF_module_mp_wsm6.F90:21