ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TI_substep_fun.H
Go to the documentation of this file.
1 #include "ERF_SrcHeaders.H"
2 
3 /**
4  * Wrapper for calling the routine that creates the fast RHS
5  */
6 auto acoustic_substepping_fun = [&](int fast_step, int n_sub, int nrk,
7  Vector<MultiFab>& S_slow_rhs,
8  const Vector<MultiFab>& S_old,
9  Vector<MultiFab>& S_stage,
10  Vector<MultiFab>& S_data,
11  const Real dtau,
12  const Real slow_dt,
13  const Real inv_fac,
14  const Real old_substep_time,
15  const Real new_substep_time)
16 {
17  BL_PROFILE("acoustic_substepping_fun");
18  if (verbose) amrex::Print() << "Fast time integration at level " << level
19  << std::setprecision(timeprecision)
20  << " from " << old_substep_time << " to " << new_substep_time
21  << " with dt = " << dtau << std::endl;
22 
23  // Define beta_s here so that it is consistent between where we make the fast coefficients
24  // and where we use them
25  // Per p2902 of Klemp-Skamarock-Dudhia-2007
26  // beta_s = -1.0 : fully explicit
27  // beta_s = 1.0 : fully implicit
28  Real beta_s;
29  if (solverChoice.substepping_type[level] == SubsteppingType::Implicit) {
30  beta_s = solverChoice.beta_s;
31  } else { // Fully explicit
32  beta_s = -1.0;
33  }
34 
35  const GpuArray<Real, AMREX_SPACEDIM> dxInv = fine_geom.InvCellSizeArray();
36 
37  // Canopy data for mom sources
38  MultiFab* forest_drag = (solverChoice.do_forest_drag) ?
39  m_forest_drag[level]->get_drag_field() : nullptr;
40 
41  // Immersed Forcing
42  MultiFab* terrain_blank = (solverChoice.terrain_type == TerrainType::ImmersedForcing ||
43  solverChoice.buildings_type == BuildingsType::ImmersedForcing ) ?
44  terrain_blanking[level].get() : nullptr;
45 
46  Vector<MultiFab> Svec_to_use;
47  if (fast_step == 0) {
48  // If this is the first substep we pass in S_old as the previous step's solution
49  Svec_to_use.push_back(MultiFab(S_old[IntVars::cons],make_alias,0,S_old[IntVars::cons].nComp()));
50  Svec_to_use.push_back(MultiFab(S_old[IntVars::xmom],make_alias,0,1));
51  Svec_to_use.push_back(MultiFab(S_old[IntVars::ymom],make_alias,0,1));
52  Svec_to_use.push_back(MultiFab(S_old[IntVars::zmom],make_alias,0,1));
53  } else {
54  // If this is not the first substep we pass in S_data as the previous step's solution
55  Svec_to_use.push_back(MultiFab(S_data[IntVars::cons],make_alias,0,S_data[IntVars::cons].nComp()));
56  Svec_to_use.push_back(MultiFab(S_data[IntVars::xmom],make_alias,0,1));
57  Svec_to_use.push_back(MultiFab(S_data[IntVars::ymom],make_alias,0,1));
58  Svec_to_use.push_back(MultiFab(S_data[IntVars::zmom],make_alias,0,1));
59  }
60 
61  make_sources(level, nrk, dtau, old_substep_time,
62  Svec_to_use, S_prim, cc_src, base_state[level], z_phys_cc[level].get(),
63  xvel_new, yvel_new,
64  qheating_rates[level].get(),
65  terrain_blank, fine_geom, solverChoice,
66  mapfac[level],
67  rhotheta_src[level].get(), rhoqt_src[level].get(),
68  dptr_wbar_sub, d_rayleigh_ptrs_at_lev,
69  d_sinesq_at_lev,
70  (solverChoice.hindcast_surface_bcs? &surface_state_interp[level] : nullptr),
71  input_sounding_data, turbPert, false);
72 
73  // *************************************************************************
74  // Set up flux registers if using two_way coupling
75  // *************************************************************************
76  const bool l_reflux = ( (solverChoice.coupling_type == CouplingType::TwoWay) && (nrk == 2) && (finest_level > 0) );
77 
78  YAFluxRegister* fr_as_crse = nullptr;
79  YAFluxRegister* fr_as_fine = nullptr;
80  if (l_reflux) {
81  if (level < finest_level) {
82  fr_as_crse = getAdvFluxReg(level+1);
83  }
84  if (level > 0) {
85  fr_as_fine = getAdvFluxReg(level);
86  }
87  }
88 
89  MultiFab base_to_use;
90  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh ) {
91  base_to_use = MultiFab(base_state_new[level], make_alias, 0, BaseState::num_comps);
92  } else {
93  base_to_use = MultiFab(base_state[level], make_alias, 0, BaseState::num_comps);
94  }
95 
96  make_mom_sources(old_substep_time, dtau, Svec_to_use,
97  z_phys_nd[level].get(), z_phys_cc[level].get(), stretched_dz_h[level],
98  xvel_new, yvel_new, zvel_new,
99  xmom_src, ymom_src, zmom_src,
100  base_to_use, forest_drag, terrain_blank,
101  cosPhi_m[level].get(), sinPhi_m[level].get(), fine_geom, solverChoice,
102  mapfac[level],
103  (solverChoice.have_geo_wind_profile) ? d_u_geos[level].data() : nullptr,
104  (solverChoice.have_geo_wind_profile) ? d_v_geos[level].data() : nullptr,
105  dptr_wbar_sub, d_rayleigh_ptrs_at_lev,
106  d_sinesq_at_lev, d_sinesq_stag_at_lev, d_sponge_ptrs_at_lev,
107  (solverChoice.hindcast_lateral_forcing? &forecast_state_interp[level] : nullptr),
108  (solverChoice.hindcast_surface_bcs? &surface_state_interp[level] : nullptr),
109  input_sounding_data, false);
110 
111  // Moving terrain
112  std::unique_ptr<MultiFab> z_t_pert;
113  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
114  {
115  z_t_pert = std::make_unique<MultiFab>(S_data[IntVars::zmom].boxArray(), S_data[IntVars::zmom].DistributionMap(), 1, 1);
116  update_terrain_substep(level, old_substep_time, new_substep_time, dtau, S_data, z_t_pert);
117  }
118 
119  bool l_use_moisture = ( solverChoice.moisture_type != MoistureType::None );
120 
121  if ( (fast_step == 0) || (solverChoice.terrain_type == TerrainType::MovingFittedMesh) )
122  {
123  make_fast_coeffs(level, fast_coeffs, S_stage, S_prim, pi_stage, fine_geom,
124  l_use_moisture, SolverChoice::mesh_type, solverChoice.gravity, solverChoice.c_p,
125  detJ_cc[level], r0, pi0, dtau, beta_s, phys_bc_type);
126  }
127 
128  bool l_rayleigh_implicit = (solverChoice.dampingChoice.rayleigh_damping_type == RayleighDampingType::FastImplicit);
129  Real l_damp_coef = solverChoice.dampingChoice.rayleigh_dampcoef;
130 
131  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
132  {
133  erf_substep_MT(fast_step, nrk, level, finest_level,
134  S_slow_rhs, Svec_to_use, S_stage, S_prim, qt, pi_stage, fast_coeffs,
135  S_data, lagged_delta_rt[level], avg_xmom[level], avg_ymom[level], avg_zmom[level],
136  cc_src, xmom_src, ymom_src, zmom_src,
137  fine_geom,
138  solverChoice.gravity, solverChoice.use_lagged_delta_rt,
139  z_t_rk[level], z_t_pert.get(),
140  z_phys_nd[level], z_phys_nd_new[level], z_phys_nd_src[level],
141  detJ_cc[level], detJ_cc_new[level], detJ_cc_src[level],
142  dtau, beta_s, inv_fac, mapfac[level],
143  fr_as_crse, fr_as_fine, l_use_moisture, l_reflux);
144 
145  } else if (solverChoice.mesh_type == MeshType::VariableDz) {
146 
147  erf_substep_T(fast_step, nrk, level, finest_level,
148  S_slow_rhs, Svec_to_use, S_stage, S_prim, qt, pi_stage, fast_coeffs,
149  S_data, lagged_delta_rt[level], avg_xmom[level], avg_ymom[level], avg_zmom[level],
150  cc_src, xmom_src, ymom_src, zmom_src,
151  fine_geom, solverChoice.gravity,
152  z_phys_nd[level], detJ_cc[level], dtau, beta_s, inv_fac,
153  mapfac[level],
154  fr_as_crse, fr_as_fine, l_use_moisture, l_reflux,
155  (l_rayleigh_implicit) ? d_sinesq_ptrs[level].data() : nullptr, l_damp_coef);
156 
157  } else {
158 
159  erf_substep_NS(fast_step, nrk, level, finest_level,
160  S_slow_rhs, Svec_to_use, S_stage, S_prim, qt, pi_stage, fast_coeffs,
161  S_data, lagged_delta_rt[level], avg_xmom[level], avg_ymom[level], avg_zmom[level],
162  cc_src, xmom_src, ymom_src, zmom_src,
163  fine_geom, solverChoice.gravity, stretched_dz_d[level],
164  dtau, beta_s, inv_fac, mapfac[level],
165  fr_as_crse, fr_as_fine, l_use_moisture, l_reflux,
166  (l_rayleigh_implicit) ? d_sinesq_ptrs[level].data() : nullptr, l_damp_coef);
167  }
168 
169  // Even if we update all the conserved variables we don't need
170  // to fillpatch the slow ones every acoustic substep
171 
172  if ( (solverChoice.vert_implicit_fac[nrk] > 0.) && !solverChoice.implicit_before_substep &&
173  (fast_step == n_sub-1) )
174  {
175  MultiFab scratch(S_data[IntVars::cons], make_alias, 0, 2);
176  MultiFab scratch_xmom(S_data[IntVars::xmom], make_alias, 0, 1);
177  MultiFab scratch_ymom(S_data[IntVars::ymom], make_alias, 0, 1);
178 #ifdef ERF_IMPLICIT_W
179  MultiFab scratch_zmom(S_data[IntVars::zmom], make_alias, 0, 1);
180 #endif
181 #include "ERF_Implicit.H"
182  }
183 
184  // Apply vertical velocity damping at the end of the timestep
185  if ( (solverChoice.dampingChoice.w_damping) &&
186  (nrk == 2) &&
187  (fast_step == n_sub-1))
188  {
189  amrex::Real w_damping_const = solverChoice.dampingChoice.w_damping_const;
190  amrex::Real w_damping_coeff = solverChoice.dampingChoice.w_damping_coeff;
191  amrex::Real cflw_lim = solverChoice.dampingChoice.w_damping_cfl;
192 
193  for ( MFIter mfi(S_data[IntVars::cons]); mfi.isValid(); ++mfi)
194  {
195  Box tbz = mfi.nodaltilebox(2);
196 
197  const Array4<const Real>& z_nd_arr = z_phys_nd[level]->const_array(mfi);
198 
199  const Array4<const Real>& cell_data = S_data[IntVars::cons].const_array(mfi);
200  const Array4< Real>& rho_w = S_data[IntVars::zmom].array(mfi);
201 
202  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
203  {
204  Real dzInv;
205  if (z_nd_arr) {
206  // average z_nd to face then diff
207  if (k == domain.smallEnd(2)) {
208  dzInv = 4.0 / (
209  z_nd_arr(i ,j ,k+1) - z_nd_arr(i ,j ,k)
210  + z_nd_arr(i+1,j ,k+1) - z_nd_arr(i+1,j ,k)
211  + z_nd_arr(i ,j+1,k+1) - z_nd_arr(i ,j+1,k)
212  + z_nd_arr(i+1,j+1,k+1) - z_nd_arr(i+1,j+1,k) );
213  } else if (k == domain.bigEnd(2)+1) {
214  dzInv = 4.0 / (
215  z_nd_arr(i ,j ,k) - z_nd_arr(i ,j ,k-1)
216  + z_nd_arr(i+1,j ,k) - z_nd_arr(i+1,j ,k-1)
217  + z_nd_arr(i ,j+1,k) - z_nd_arr(i ,j+1,k-1)
218  + z_nd_arr(i+1,j+1,k) - z_nd_arr(i+1,j+1,k-1) );
219  } else {
220  // dz = 0.5 * (dz[i,j,k] + dz[i,j,k+1])
221  // = 0.5 * (z_nd_face[i,j,k+1] - z_nd_face[i,j,k-1])
222  dzInv = 8.0 / (
223  z_nd_arr(i ,j ,k+1) - z_nd_arr(i ,j ,k-1)
224  + z_nd_arr(i+1,j ,k+1) - z_nd_arr(i+1,j ,k-1)
225  + z_nd_arr(i ,j+1,k+1) - z_nd_arr(i ,j+1,k-1)
226  + z_nd_arr(i+1,j+1,k+1) - z_nd_arr(i+1,j+1,k-1) );
227  }
228  } else {
229  dzInv = dxInv[2];
230  }
231 
232  Real dampcoef = (w_damping_const > 0)
233  ? w_damping_const
234  : w_damping_coeff / (dzInv * dt_advance * dt_advance);
235 
236  Real rho_on_w_face = 0.5 * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp) );
237  Real wmag = std::abs(rho_w(i,j,k) / rho_on_w_face);
238 
239  // This is the dynamics timestep, not the slow or fast dt
240  Real cflw = wmag * dt_advance * dzInv;
241 
242  if (cflw > cflw_lim) {
243  Real sgn_w = (rho_w(i,j,k) > 0) ? 1.0 : -1.0;
244  rho_w(i, j, k) -= rho_on_w_face * sgn_w * dampcoef * (cflw - cflw_lim) * dt_advance;
245  // Note: These prints are not immediately flushed to screen
246 #ifdef AMREX_USE_GPU
247  AMREX_DEVICE_PRINTF("t=%f w-damping applied at (%d,%d,%d) for w-CFL = %f > %f : %f --> %f\n",
248  old_time+dt_advance, i,j,k, cflw, cflw_lim, sgn_w*wmag, rho_w(i,j,k)/rho_on_w_face);
249 #else
250  printf("t=%f w-damping applied at (%d,%d,%d) for w-CFL = %f > %f : %f --> %f\n",
251  old_time+dt_advance, i,j,k, cflw, cflw_lim, sgn_w*wmag, rho_w(i,j,k)/rho_on_w_face);
252 #endif
253  }
254  });
255  }
256  }
257 
258  // NOTE: Numerical diffusion is tested on in FillPatchIntermediate and dictates the size of the
259  // box over which VelocityToMomentum is computed. V2M requires one more ghost cell be
260  // filled for rho than velocity. This logical condition ensures we fill enough ghost cells
261  // when use_num_diff is true.
262  int ng_cons = S_data[IntVars::cons].nGrowVect().max() - 1;
263  int ng_vel = S_data[IntVars::xmom].nGrowVect().max();
264  if (!solverChoice.use_num_diff) {
265  ng_cons = 1;
266  ng_vel = 1;
267  }
268  apply_bcs(S_data, new_substep_time, ng_cons, ng_vel, fast_only=true, vel_and_mom_synced=false);
269 };
if(l_use_mynn &&start_comp<=RhoKE_comp &&end_comp >=RhoKE_comp)
Definition: ERF_AddQKESources.H:2
#define Rho_comp
Definition: ERF_IndexDefines.H:36
void make_fast_coeffs(int, MultiFab &fast_coeffs, Vector< MultiFab > &S_stage_data, const MultiFab &S_stage_prim, const MultiFab &pi_stage, const amrex::Geometry geom, bool l_use_moisture, MeshType mesh_type, Real gravity, Real c_p, std::unique_ptr< MultiFab > &detJ_cc, const MultiFab *r0, const MultiFab *pi0, Real dtau, Real beta_s, amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
Definition: ERF_MakeFastCoeffs.cpp:26
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, 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:33
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void erf_substep_MT(int step, int, int level, int finest_level, Vector< MultiFab > &S_slow_rhs, const Vector< MultiFab > &S_prev, Vector< MultiFab > &S_stg_data, const MultiFab &S_stg_prim, const MultiFab &qt, const MultiFab &pi_stage, const MultiFab &fast_coeffs, Vector< MultiFab > &S_data, MultiFab &lagged_delta_rt, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &cc_src, const MultiFab &xmom_src, const MultiFab &ymom_src, const MultiFab &zmom_src, const Geometry geom, const Real gravity, const bool use_lagged_delta_rt, std::unique_ptr< MultiFab > &z_t_rk, const MultiFab *z_t_pert, std::unique_ptr< MultiFab > &z_phys_nd_old, std::unique_ptr< MultiFab > &z_phys_nd_new, std::unique_ptr< MultiFab > &z_phys_nd_stg, std::unique_ptr< MultiFab > &detJ_cc_old, std::unique_ptr< MultiFab > &detJ_cc_new, std::unique_ptr< MultiFab > &detJ_cc_stg, const Real dtau, const Real beta_s, const Real facinv, Vector< std::unique_ptr< MultiFab >> &mapfac, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine, bool l_use_moisture, bool l_reflux)
Definition: ERF_Substep_MT.cpp:49
void erf_substep_NS(int step, int nrk, int level, int finest_level, Vector< MultiFab > &S_slow_rhs, const Vector< MultiFab > &S_prev, Vector< MultiFab > &S_stage_data, const MultiFab &S_stage_prim, const MultiFab &qt, const MultiFab &pi_stage, const MultiFab &fast_coeffs, Vector< MultiFab > &S_data, MultiFab &lagged_delta_rt, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &cc_src, const MultiFab &xmom_src, const MultiFab &ymom_src, const MultiFab &zmom_src, const Geometry geom, const Real gravity, amrex::Gpu::DeviceVector< amrex::Real > &stretched_dz_d, const Real dtau, const Real beta_s, const Real facinv, Vector< std::unique_ptr< MultiFab >> &mapfac, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine, bool l_use_moisture, bool l_reflux, const amrex::Real *sinesq_stag_d, const Real l_damp_coef)
Definition: ERF_Substep_NS.cpp:42
void erf_substep_T(int step, int, int level, int finest_level, Vector< MultiFab > &S_slow_rhs, const Vector< MultiFab > &S_prev, Vector< MultiFab > &S_stage_data, const MultiFab &S_stage_prim, const MultiFab &qt, const MultiFab &pi_stage, const MultiFab &fast_coeffs, Vector< MultiFab > &S_data, MultiFab &lagged_delta_rt, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &cc_src, const MultiFab &xmom_src, const MultiFab &ymom_src, const MultiFab &zmom_src, const Geometry geom, const Real gravity, std::unique_ptr< MultiFab > &z_phys_nd, std::unique_ptr< MultiFab > &detJ_cc, const Real dtau, const Real beta_s, const Real facinv, Vector< std::unique_ptr< MultiFab >> &mapfac, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine, bool l_use_moisture, bool l_reflux, const Real *sinesq_stag_d, const Real l_damp_coef)
Definition: ERF_Substep_T.cpp:43
auto acoustic_substepping_fun
Definition: ERF_TI_substep_fun.H:6
auto apply_bcs
Definition: ERF_TI_utils.H:73
auto update_terrain_substep
Definition: ERF_TI_utils.H:238
@ num_comps
Definition: ERF_IndexDefines.H:68
@ 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
static MeshType mesh_type
Definition: ERF_DataStruct.H:1048