ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Substep_MT.cpp File Reference
Include dependency graph for ERF_Substep_MT.cpp:

Functions

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, bool, const Real *sinesq_stag_d, const Real l_damp_coef)
 

Function Documentation

◆ erf_substep_MT()

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,
bool  ,
const Real sinesq_stag_d,
const Real  l_damp_coef 
)

Function for computing the fast RHS with moving terrain

Parameters
[in]stepwhich fast time step within each Runge-Kutta step
[in]nrkwhich Runge-Kutta step
[in]levellevel of resolution
[in]finest_levelfinest level of resolution
[in]S_slow_rhsslow RHS computed in erf_slow_rhs_pre
[in]S_prevprevious solution
[in]S_stg_datasolution at previous RK stage
[in]S_stg_primprimitive variables at previous RK stage
[in]pi_stageExner function at previous RK stage
[in]fast_coeffscoefficients for the tridiagonal solve used in the fast integrator
[out]S_datacurrent solution
[in,out]lagged_delta_rt
[in,out]avg_xmomtime-averaged x-momentum to be used for updating slow variables
[in,out]avg_ymomtime-averaged y-momentum to be used for updating slow variables
[in,out]avg_zmomtime-averaged z-momentum to be used for updating slow variables
[in]cc_srcsource terms for conserved variables
[in]xmom_srcsource terms for x-momentum
[in]ymom_srcsource terms for y-momentum
[in]zmom_srcsource terms for z-momentum
[in]geomcontainer for geometric information
[in]gravityMagnitude of gravity
[in]use_lagged_delta_rtdefine lagged_delta_rt for our next step
[in]z_t_rkrate of change of grid height – only relevant for moving terrain
[in]z_t_pertrate of change of grid height – interpolated between RK stages
[in]z_phys_nd_oldheight coordinate at nodes at old time
[in]z_phys_nd_newheight coordinate at nodes at new time
[in]z_phys_nd_stgheight coordinate at nodes at previous stage
[in]detJ_cc_oldJacobian of the metric transformation at old time
[in]detJ_cc_newJacobian of the metric transformation at new time
[in]detJ_cc_stgJacobian of the metric transformation at previous stage
[in]dtaufast time step
[in]beta_sCoefficient which determines how implicit vs explicit the solve is
[in]facinvinverse factor for time-averaging the momenta
[in]mapfacvector of map 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
[in]l_use_moisture
[in]l_refluxshould we add fluxes to the FluxRegisters?
[in]l_damp_coef
89 {
90  BL_PROFILE_REGION("erf_substep_MT()");
91 
92  Real beta_1 = myhalf * (one - beta_s); // multiplies explicit terms
93  Real beta_2 = myhalf * (one + beta_s); // multiplies implicit terms
94 
95  // How much do we project forward the (rho theta) that is used in the horizontal momentum equations
96  Real beta_d = Real(0.1);
97 
98  Real RvOverRd = R_v / R_d;
99 
100  bool l_rayleigh_impl_for_w = (sinesq_stag_d != nullptr);
101 
102  const Real* dx = geom.CellSize();
103  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
104 
105  Real dxi = dxInv[0];
106  Real dyi = dxInv[1];
107  Real dzi = dxInv[2];
108 
109  MultiFab coeff_A_mf(fast_coeffs, make_alias, 0, 1);
110  MultiFab inv_coeff_B_mf(fast_coeffs, make_alias, 1, 1);
111  MultiFab coeff_C_mf(fast_coeffs, make_alias, 2, 1);
112  MultiFab coeff_P_mf(fast_coeffs, make_alias, 3, 1);
113  MultiFab coeff_Q_mf(fast_coeffs, make_alias, 4, 1);
114 
115  // *************************************************************************
116  // Set gravity as a vector
117  const Array<Real,AMREX_SPACEDIM> grav{zero, zero, -gravity};
118  const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
119 
120  MultiFab extrap(S_data[IntVars::cons].boxArray(),S_data[IntVars::cons].DistributionMap(),1,1);
121 
122  MultiFab Omega(S_data[IntVars::zmom].boxArray(), S_data[IntVars::zmom].DistributionMap(), 1, 1);
123 
124  // *************************************************************************
125  // Define updates in the current RK stg
126  // *************************************************************************
127 #ifdef _OPENMP
128 #pragma omp parallel if (Gpu::notInLaunchRegion())
129 #endif
130  {
131  FArrayBox temp_rhs_fab;
132 
133  FArrayBox RHS_fab;
134  FArrayBox soln_fab;
135 
136  std::array<FArrayBox,AMREX_SPACEDIM> flux;
137 
138  // NOTE: we leave tiling off here for efficiency -- to make this loop work with tiling
139  // will require additional changes
140  for ( MFIter mfi(S_stg_data[IntVars::cons],false); mfi.isValid(); ++mfi)
141  {
142  Box bx = mfi.tilebox();
143  Box tbx = surroundingNodes(bx,0);
144  Box tby = surroundingNodes(bx,1);
145  Box tbz = surroundingNodes(bx,2);
146 
147  Box vbx = mfi.validbox();
148  const auto& vbx_hi = ubound(vbx);
149 
150  const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
151  const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
152  const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
153  const Array4<Real const>& cc_src_arr = cc_src.const_array(mfi);
154 
155  const Array4<const Real> & stg_cons = S_stg_data[IntVars::cons].const_array(mfi);
156  const Array4<const Real> & stg_xmom = S_stg_data[IntVars::xmom].const_array(mfi);
157  const Array4<const Real> & stg_ymom = S_stg_data[IntVars::ymom].const_array(mfi);
158  const Array4<const Real> & stg_zmom = S_stg_data[IntVars::zmom].const_array(mfi);
159  const Array4<const Real> & prim = S_stg_prim.const_array(mfi);
160  const Array4<const Real> & qt_arr = qt.const_array(mfi);
161 
162  const Array4<const Real>& slow_rhs_cons = S_slow_rhs[IntVars::cons].const_array(mfi);
163  const Array4<const Real>& slow_rhs_rho_u = S_slow_rhs[IntVars::xmom].const_array(mfi);
164  const Array4<const Real>& slow_rhs_rho_v = S_slow_rhs[IntVars::ymom].const_array(mfi);
165  const Array4<const Real>& slow_rhs_rho_w = S_slow_rhs[IntVars::zmom].const_array(mfi);
166 
167  const Array4<Real>& cur_cons = S_data[IntVars::cons].array(mfi);
168  const Array4<Real>& cur_xmom = S_data[IntVars::xmom].array(mfi);
169  const Array4<Real>& cur_ymom = S_data[IntVars::ymom].array(mfi);
170  const Array4<Real>& cur_zmom = S_data[IntVars::zmom].array(mfi);
171 
172  const Array4<Real>& lagged = lagged_delta_rt.array(mfi);
173 
174  const Array4<const Real>& prev_cons = S_prev[IntVars::cons].const_array(mfi);
175  const Array4<const Real>& prev_xmom = S_prev[IntVars::xmom].const_array(mfi);
176  const Array4<const Real>& prev_ymom = S_prev[IntVars::ymom].const_array(mfi);
177  const Array4<const Real>& prev_zmom = S_prev[IntVars::zmom].const_array(mfi);
178 
179  // These store the advection momenta which we will use to update the slow variables
180  const Array4<Real>& avg_xmom_arr = avg_xmom.array(mfi);
181  const Array4<Real>& avg_ymom_arr = avg_ymom.array(mfi);
182  const Array4<Real>& avg_zmom_arr = avg_zmom.array(mfi);
183 
184  const Array4<const Real>& z_nd_old = z_phys_nd_old->const_array(mfi);
185  const Array4<const Real>& z_nd_new = z_phys_nd_new->const_array(mfi);
186  const Array4<const Real>& z_nd_stg = z_phys_nd_stg->const_array(mfi);
187  const Array4<const Real>& detJ_old = detJ_cc_old->const_array(mfi);
188  const Array4<const Real>& detJ_new = detJ_cc_new->const_array(mfi);
189  const Array4<const Real>& detJ_stg = detJ_cc_stg->const_array(mfi);
190 
191  const Array4<const Real>& z_t_arr = z_t_rk->const_array(mfi);
192  const Array4<const Real>& zp_t_arr = z_t_pert->const_array(mfi);
193 
194  const Array4< Real>& omega_arr = Omega.array(mfi);
195 
196  // Map factors
197  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
198  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
199  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
200  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
201 
202  // *********************************************************************
203  // This must be done before we set cur_xmom and cur_ymom, since those
204  // in fact point to the same array as prev_xmom and prev_ymom
205  // *********************************************************************
206  Box gbxo = mfi.nodaltilebox(2);
207  {
208  BL_PROFILE("fast_MT_making_omega");
209  Box gbxo_lo = gbxo; gbxo_lo.setBig(2,0);
210  ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
211  omega_arr(i,j,k) = zero;
212  });
213  Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
214  ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
215  omega_arr(i,j,k) = prev_zmom(i,j,k) - stg_zmom(i,j,k) - zp_t_arr(i,j,k);
216  });
217  Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
218  ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
219  omega_arr(i,j,k) =
220  ( OmegaFromW(i,j,k,prev_zmom(i,j,k),prev_xmom,prev_ymom,mf_ux,mf_vy,z_nd_old,dxInv)
221  -OmegaFromW(i,j,k, stg_zmom(i,j,k), stg_xmom, stg_ymom,mf_ux,mf_vy,z_nd_old,dxInv) )
222  - zp_t_arr(i,j,k);
223  });
224  } // end profile
225  // *********************************************************************
226 
227  const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
228 
229  const Array4<Real>& theta_extrap = extrap.array(mfi);
230 
231  // Note: it is important to grow the tilebox rather than use growntilebox because
232  // we need to fill the ghost cells of the tilebox so we can use them below
233  Box gbx = mfi.tilebox(); gbx.grow(1);
234  Box gtbx = mfi.nodaltilebox(0); gtbx.grow(1); gtbx.setSmall(2,0);
235  Box gtby = mfi.nodaltilebox(1); gtby.grow(1); gtby.setSmall(2,0);
236 
237  if (step == 0) {
238  ParallelFor(gbx,
239  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
240  cur_cons(i,j,k,Rho_comp) = prev_cons(i,j,k,Rho_comp);
241  cur_cons(i,j,k,RhoTheta_comp) = prev_cons(i,j,k,RhoTheta_comp);
242 
243  Real delta_rt = cur_cons(i,j,k,RhoTheta_comp) - stg_cons(i,j,k,RhoTheta_comp);
244  theta_extrap(i,j,k) = delta_rt;
245 
246  // NOTE: qv is not changing over the fast steps so we use the stage data
247  Real qv = (l_use_moisture) ? prim(i,j,k,PrimQ1_comp) : zero;
248  theta_extrap(i,j,k) *= (one + RvOverRd*qv);
249 
250  // We define lagged_delta_rt for our next step as the current delta_rt
251  lagged(i,j,k) = delta_rt;
252  });
253  } else if (use_lagged_delta_rt) {
254  // This is the default for cases with no or static terrain
255  ParallelFor(gbx,
256  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
257  Real delta_rt = cur_cons(i,j,k,RhoTheta_comp) - stg_cons(i,j,k,RhoTheta_comp);
258  theta_extrap(i,j,k) = delta_rt + beta_d * (delta_rt - lagged(i,j,k));
259 
260  // NOTE: qv is not changing over the fast steps so we use the stage data
261  Real qv = (l_use_moisture) ? prim(i,j,k,PrimQ1_comp) : zero;
262  theta_extrap(i,j,k) *= (one + RvOverRd*qv);
263 
264  // We define lagged_delta_rt for our next step as the current delta_rt
265  lagged(i,j,k) = delta_rt;
266  });
267  } else {
268  // For the moving wave problem, this choice seems more robust
269  ParallelFor(gbx,
270  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
271  theta_extrap(i,j,k) = cur_cons(i,j,k,RhoTheta_comp) - stg_cons(i,j,k,RhoTheta_comp);
272 
273  // NOTE: qv is not changing over the fast steps so we use the stage data
274  Real qv = (l_use_moisture) ? prim(i,j,k,PrimQ1_comp) : zero;
275  theta_extrap(i,j,k) *= (one + RvOverRd*qv);
276  });
277  } // if step
278 
279  RHS_fab.resize (tbz,1, The_Async_Arena());
280  soln_fab.resize (tbz,1, The_Async_Arena());
281  temp_rhs_fab.resize(tbz,2, The_Async_Arena());
282 
283  auto const& RHS_a = RHS_fab.array();
284  auto const& soln_a = soln_fab.array();
285  auto const& temp_rhs_arr = temp_rhs_fab.array();
286 
287  auto const& coeffA_a = coeff_A_mf.array(mfi);
288  auto const& inv_coeffB_a = inv_coeff_B_mf.array(mfi);
289  auto const& coeffC_a = coeff_C_mf.array(mfi);
290  auto const& coeffP_a = coeff_P_mf.array(mfi);
291  auto const& coeffQ_a = coeff_Q_mf.array(mfi);
292 
293  // *********************************************************************
294  // Define updates in the RHS of {x, y, z}-momentum equations
295  // *********************************************************************
296  {
297  BL_PROFILE("substep_xymom_T");
298  ParallelFor(tbx, tby,
299  [=] AMREX_GPU_DEVICE (int i, int j, int k)
300  {
301  // Add (negative) gradient of (rho theta) multiplied by lagged "pi"
302  Real h_xi_old = Compute_h_xi_AtIface(i, j, k, dxInv, z_nd_old);
303  Real h_zeta_old = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd_old);
304  Real gp_xi = (theta_extrap(i,j,k) - theta_extrap(i-1,j,k)) * dxi;
305  Real gp_zeta_on_iface = (k == 0) ?
306  myhalf * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
307  - theta_extrap(i-1,j,k ) - theta_extrap(i,j,k ) ) :
308  fourth * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
309  - theta_extrap(i-1,j,k-1) - theta_extrap(i,j,k-1) );
310  Real gpx = h_zeta_old * gp_xi - h_xi_old * gp_zeta_on_iface;
311  gpx *= mf_ux(i,j,0);
312 
313  Real q = (l_use_moisture) ? myhalf * (qt_arr(i-1,j,k) + qt_arr(i,j,k)) : zero;
314 
315  Real pi_c = myhalf * (pi_stage_ca(i-1,j,k) + pi_stage_ca(i ,j,k));
316  Real fast_rhs_rho_u = -Gamma * R_d * pi_c * gpx / (one + q);
317 
318  // We have already scaled the source terms to have the extra factor of dJ
319  cur_xmom(i,j,k) = h_zeta_old * prev_xmom(i,j,k) + dtau * fast_rhs_rho_u
320  + dtau * slow_rhs_rho_u(i,j,k)
321  + dtau * xmom_src_arr(i,j,k);
322  },
323  [=] AMREX_GPU_DEVICE (int i, int j, int k)
324  {
325  // Add (negative) gradient of (rho theta) multiplied by lagged "pi"
326  Real h_eta_old = Compute_h_eta_AtJface(i, j, k, dxInv, z_nd_old);
327  Real h_zeta_old = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd_old);
328  Real gp_eta = (theta_extrap(i,j,k) -theta_extrap(i,j-1,k)) * dyi;
329  Real gp_zeta_on_jface = (k == 0) ?
330  myhalf * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
331  - theta_extrap(i,j,k ) - theta_extrap(i,j-1,k ) ) :
332  fourth * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
333  - theta_extrap(i,j,k-1) - theta_extrap(i,j-1,k-1) );
334  Real gpy = h_zeta_old * gp_eta - h_eta_old * gp_zeta_on_jface;
335  gpy *= mf_vy(i,j,0);
336 
337  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j-1,k) + qt_arr(i,j,k)) : zero;
338 
339  Real pi_c = myhalf * (pi_stage_ca(i,j-1,k) + pi_stage_ca(i,j ,k));
340  Real fast_rhs_rho_v = -Gamma * R_d * pi_c * gpy / (one + q);
341 
342  // We have already scaled the source terms to have the extra factor of dJ
343  cur_ymom(i, j, k) = h_zeta_old * prev_ymom(i,j,k) + dtau * fast_rhs_rho_v
344  + dtau * slow_rhs_rho_v(i,j,k)
345  + dtau * ymom_src_arr(i,j,k);
346  });
347  } // end profile
348 
349  // *************************************************************************
350  // Define flux arrays for use in advection
351  // *************************************************************************
352  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
353  flux[dir].resize(surroundingNodes(bx,dir),2,The_Async_Arena());
354  flux[dir].setVal<RunOn::Device>(0);
355  }
356  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
357  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
358 
359  // *********************************************************************
360  {
361  BL_PROFILE("fast_T_making_rho_rhs");
362  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
363  {
364  Real h_zeta_stg_xlo = Compute_h_zeta_AtIface(i, j , k, dxInv, z_nd_stg);
365  Real h_zeta_stg_xhi = Compute_h_zeta_AtIface(i+1,j , k, dxInv, z_nd_stg);
366  Real xflux_lo = cur_xmom(i ,j,k) - stg_xmom(i ,j,k)*h_zeta_stg_xlo;
367  Real xflux_hi = cur_xmom(i+1,j,k) - stg_xmom(i+1,j,k)*h_zeta_stg_xhi;
368 
369  Real h_zeta_stg_yhi = Compute_h_zeta_AtJface(i, j+1, k, dxInv, z_nd_stg);
370  Real h_zeta_stg_ylo = Compute_h_zeta_AtJface(i, j , k, dxInv, z_nd_stg);
371  Real yflux_lo = cur_ymom(i,j ,k) - stg_ymom(i,j ,k)*h_zeta_stg_ylo;
372  Real yflux_hi = cur_ymom(i,j+1,k) - stg_ymom(i,j+1,k)*h_zeta_stg_yhi;
373 
374  // NOTE: we are saving the (1/J) weighting for later when we add this to rho and theta
375  temp_rhs_arr(i,j,k,0) = ( xflux_hi - xflux_lo ) * dxi + ( yflux_hi - yflux_lo ) * dyi;
376  temp_rhs_arr(i,j,k,1) = (( xflux_hi * (prim(i,j,k,0) + prim(i+1,j,k,0)) -
377  xflux_lo * (prim(i,j,k,0) + prim(i-1,j,k,0)) ) * dxi +
378  ( yflux_hi * (prim(i,j,k,0) + prim(i,j+1,k,0)) -
379  yflux_lo * (prim(i,j,k,0) + prim(i,j-1,k,0)) ) * dyi) * myhalf;
380 
381  if (l_reflux) {
382  (flx_arr[0])(i,j,k,0) = xflux_lo;
383  (flx_arr[0])(i,j,k,1) = (flx_arr[0])(i ,j,k,0) * myhalf * (prim(i,j,k,0) + prim(i-1,j,k,0));
384 
385  (flx_arr[1])(i,j,k,0) = yflux_lo;
386  (flx_arr[1])(i,j,k,1) = (flx_arr[1])(i,j ,k,0) * myhalf * (prim(i,j,k,0) + prim(i,j-1,k,0));
387 
388  if (i == vbx_hi.x) {
389  (flx_arr[0])(i+1,j,k,0) = xflux_hi;
390  (flx_arr[0])(i+1,j,k,1) = (flx_arr[0])(i+1,j,k,0) * myhalf * (prim(i,j,k,0) + prim(i+1,j,k,0));
391  }
392  if (j == vbx_hi.y) {
393  (flx_arr[1])(i,j+1,k,0) = yflux_hi;
394  (flx_arr[1])(i,j+1,k,1) = (flx_arr[1])(i,j+1,k,0) * myhalf * (prim(i,j,k,0) + prim(i,j+1,k,0));
395  }
396  }
397  });
398  } // end profile
399 
400  ParallelFor(tbx, tby,
401  [=] AMREX_GPU_DEVICE (int i, int j, int k)
402  {
403  Real h_zeta_new = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd_new);
404  cur_xmom(i, j, k) /= h_zeta_new;
405  avg_xmom_arr(i,j,k) += facinv*(cur_xmom(i,j,k) - stg_xmom(i,j,k));
406  },
407  [=] AMREX_GPU_DEVICE (int i, int j, int k)
408  {
409  Real h_zeta_new = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd_new);
410  cur_ymom(i, j, k) /= h_zeta_new;
411  avg_ymom_arr(i,j,k) += facinv*(cur_ymom(i,j,k) - stg_ymom(i,j,k));
412  });
413 
414  Box bx_shrunk_in_k = bx;
415  int klo = tbz.smallEnd(2);
416  int khi = tbz.bigEnd(2);
417  bx_shrunk_in_k.setSmall(2,klo+1);
418  bx_shrunk_in_k.setBig(2,khi-1);
419 
420  // Note that the notes use "g" to mean the magnitude of gravity, so it is positive
421  // We set grav_gpu[2] to be the vector component which is negative
422  // We define halfg to match the notes (which is why we take the absolute value)
423  Real halfg = std::abs(myhalf * grav_gpu[2]);
424 
425  {
426  BL_PROFILE("fast_loop_on_shrunk_t");
427  //Note we don't act on the bottom or top boundaries of the domain
428  ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (int i, int j, int k)
429  {
430  Real dJ_old_kface = myhalf * (detJ_old(i,j,k) + detJ_old(i,j,k-1));
431  Real dJ_new_kface = myhalf * (detJ_new(i,j,k) + detJ_new(i,j,k-1));
432  Real dJ_stg_kface = myhalf * (detJ_stg(i,j,k) + detJ_stg(i,j,k-1));
433 
434  Real coeff_P = coeffP_a(i,j,k);
435  Real coeff_Q = coeffQ_a(i,j,k);
436 
437  Real theta_t_lo = myhalf * ( prim(i,j,k-2,PrimTheta_comp) + prim(i,j,k-1,PrimTheta_comp) );
438  Real theta_t_mid = myhalf * ( prim(i,j,k-1,PrimTheta_comp) + prim(i,j,k ,PrimTheta_comp) );
439  Real theta_t_hi = myhalf * ( prim(i,j,k ,PrimTheta_comp) + prim(i,j,k+1,PrimTheta_comp) );
440 
441  // line 2 last two terms (order dtau)
442  Real R0_tmp = coeff_P * cur_cons(i,j,k ,RhoTheta_comp)
443  + coeff_Q * cur_cons(i,j,k-1,RhoTheta_comp)
444  - coeff_P * stg_cons(i,j,k ,RhoTheta_comp) * (dJ_stg_kface/dJ_old_kface)
445  - coeff_Q * stg_cons(i,j,k-1,RhoTheta_comp) * (dJ_stg_kface/dJ_old_kface)
446  - halfg * ( cur_cons(i,j,k,Rho_comp) + cur_cons(i,j,k-1,Rho_comp) )
447  + halfg * ( stg_cons(i,j,k,Rho_comp) + stg_cons(i,j,k-1,Rho_comp) ) * (dJ_stg_kface/dJ_old_kface);
448 
449  // line 3 residuals (order dtau^2) one <-> beta_2
450  Real R1_tmp = - halfg * ( slow_rhs_cons(i,j,k ,Rho_comp) + slow_rhs_cons(i,j,k-1,Rho_comp))
451  + coeff_P * slow_rhs_cons(i,j,k ,RhoTheta_comp) + coeff_Q * slow_rhs_cons(i,j,k-1,RhoTheta_comp);
452 
453  Real Omega_kp1 = omega_arr(i,j,k+1);
454  Real Omega_k = omega_arr(i,j,k );
455  Real Omega_km1 = omega_arr(i,j,k-1);
456 
457  Real detJdiff = (detJ_old(i,j,k) - detJ_old(i,j,k-1)) / (detJ_old(i,j,k)*detJ_old(i,j,k-1));
458 
459  // consolidate lines 4&5 (order dtau^2)
460  R1_tmp += halfg * ( beta_1 * dzi * (Omega_kp1/detJ_old(i,j,k) + detJdiff*Omega_k - Omega_km1/detJ_old(i,j,k-1))
461  + temp_rhs_arr(i,j,k,Rho_comp)/detJ_old(i,j,k) + temp_rhs_arr(i,j,k-1,Rho_comp)/detJ_old(i,j,k-1) );
462 
463  // consolidate lines 6&7 (order dtau^2)
464  R1_tmp += -(
465  coeff_P/detJ_old(i,j,k ) * ( beta_1 * dzi * (Omega_kp1*theta_t_hi - Omega_k*theta_t_mid)
466  +temp_rhs_arr(i,j,k ,RhoTheta_comp) ) +
467  coeff_Q/detJ_old(i,j,k-1) * ( beta_1 * dzi * (Omega_k*theta_t_mid - Omega_km1*theta_t_lo)
468  +temp_rhs_arr(i,j,k-1,RhoTheta_comp) ) );
469 
470  // line 1
471  RHS_a(i,j,k) = prev_zmom(i,j,k) - (dJ_stg_kface/dJ_old_kface) * stg_zmom(i,j,k)
472  + dtau * slow_rhs_rho_w(i,j,k) / dJ_stg_kface
473  + dtau * zmom_src_arr(i,j,k);
474 
475  RHS_a(i,j,k) += dtau * R0_tmp;
476 
477  RHS_a(i,j,k) += dtau * dtau*beta_2*R1_tmp;
478 
479  // We cannot use omega_arr here since that was built with old_rho_u and old_rho_v ...
480  Real UppVpp = (dJ_new_kface/dJ_old_kface) * OmegaFromW(i,j,k,0.,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
481  -(dJ_stg_kface/dJ_old_kface) * OmegaFromW(i,j,k,0.,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
482  RHS_a(i,j,k) += UppVpp;
483  });
484  } // end profile
485 
486  Box b2d = tbz; // Copy constructor
487  b2d.setRange(2,0);
488 
489  auto const lo = lbound(bx);
490  auto const hi = ubound(bx);
491 
492  {
493  BL_PROFILE("substep_b2d_loop_t");
494 
495 #ifdef AMREX_USE_GPU
496  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int)
497  {
498  // Moving terrain
499  Real rho_on_bdy = myhalf * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
500  RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,lo.z);
501 
502  soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
503 
504  RHS_a(i,j,hi.z+1) = dtau * (slow_rhs_rho_w(i,j,hi.z+1) + zmom_src_arr(i,j,hi.z+1));
505 
506  for (int k = lo.z+1; k <= hi.z+1; k++) {
507  soln_a(i,j,k) = (RHS_a(i,j,k)-coeffA_a(i,j,k)*soln_a(i,j,k-1)) * inv_coeffB_a(i,j,k);
508  }
509 
510  for (int k = hi.z; k >= lo.z; k--) {
511  soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
512  }
513 
514  // We assume that Omega == w at the top boundary and that changes in J there are irrelevant
515  cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
516  });
517 #else
518  for (int j = lo.y; j <= hi.y; ++j) {
519  AMREX_PRAGMA_SIMD
520  for (int i = lo.x; i <= hi.x; ++i) {
521 
522  Real rho_on_bdy = myhalf * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
523  RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,lo.z);
524 
525  soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
526  }
527  }
528 
529  for (int j = lo.y; j <= hi.y; ++j) {
530  AMREX_PRAGMA_SIMD
531  for (int i = lo.x; i <= hi.x; ++i) {
532  RHS_a(i,j,hi.z+1) = dtau * (slow_rhs_rho_w(i,j,hi.z+1) + zmom_src_arr(i,j,hi.z+1));
533  }
534  }
535  for (int k = lo.z+1; k <= hi.z+1; ++k) {
536  for (int j = lo.y; j <= hi.y; ++j) {
537  AMREX_PRAGMA_SIMD
538  for (int i = lo.x; i <= hi.x; ++i) {
539  soln_a(i,j,k) = (RHS_a(i,j,k)-coeffA_a(i,j,k)*soln_a(i,j,k-1)) * inv_coeffB_a(i,j,k);
540  }
541  }
542  }
543  for (int k = hi.z; k >= lo.z; --k) {
544  for (int j = lo.y; j <= hi.y; ++j) {
545  AMREX_PRAGMA_SIMD
546  for (int i = lo.x; i <= hi.x; ++i) {
547  soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
548  }
549  }
550  }
551 
552  // We assume that Omega == w at the top boundary and that changes in J there are irrelevant
553  for (int j = lo.y; j <= hi.y; ++j) {
554  AMREX_PRAGMA_SIMD
555  for (int i = lo.x; i <= hi.x; ++i) {
556  cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
557  }
558  }
559 #endif
560  } // end profile
561 
562  {
563  BL_PROFILE("substep_new_drhow");
564  tbz.setBig(2,hi.z);
565  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
566  {
567  Real rho_on_face = myhalf * (cur_cons(i,j,k,Rho_comp) + cur_cons(i,j,k-1,Rho_comp));
568 
569  if (k == lo.z) {
570  cur_zmom(i,j,k) = WFromOmega(i,j,k,rho_on_face*(z_t_arr(i,j,k)+zp_t_arr(i,j,k)),
571  cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv);
572 
573  // We need to set this here because it is used to define zflux_lo below
574  soln_a(i,j,k) = zero;
575 
576  } else {
577 
578  Real UppVpp = WFromOmega(i,j,k,zero,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
579  - WFromOmega(i,j,k,zero,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
580  Real wpp = soln_a(i,j,k) + UppVpp;
581  Real dJ_old_kface = myhalf * (detJ_old(i,j,k) + detJ_old(i,j,k-1));
582  Real dJ_new_kface = myhalf * (detJ_new(i,j,k) + detJ_new(i,j,k-1));
583 
584  cur_zmom(i,j,k) = dJ_old_kface * (stg_zmom(i,j,k) + wpp);
585  cur_zmom(i,j,k) /= dJ_new_kface;
586 
587  soln_a(i,j,k) = OmegaFromW(i,j,k,cur_zmom(i,j,k),cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
588  - OmegaFromW(i,j,k,stg_zmom(i,j,k),stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
589  soln_a(i,j,k) -= rho_on_face * zp_t_arr(i,j,k);
590  }
591 
592  if (l_rayleigh_impl_for_w && k > 0) {
593  Real damping_coeff = l_damp_coef * dtau * sinesq_stag_d[k];
594  cur_zmom(i,j,k) /= (one + damping_coeff);
595  }
596  });
597  } // end profile
598 
599  // **************************************************************************
600  // Define updates in the RHS of rho and (rho theta)
601  // **************************************************************************
602  {
603  BL_PROFILE("fast_rho_final_update");
604  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
605  {
606  Real zflux_lo = beta_2 * soln_a(i,j,k ) + beta_1 * omega_arr(i,j,k);
607  Real zflux_hi = beta_2 * soln_a(i,j,k+1) + beta_1 * omega_arr(i,j,k+1);
608 
609  // Note that in the solve we effectively impose new_drho_w(i,j,vbx_hi.z+1)=0
610  // so we don't update avg_zmom at k=vbx_hi.z+1
611  avg_zmom_arr(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
612  if (l_reflux) {
613  (flx_arr[2])(i,j,k,0) = zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
614  }
615 
616  // Note that the factor of (1/J) in the fast source term is canceled
617  // when we multiply old and new by detJ_old and detJ_new , respectively
618  // We have already scaled the slow source term to have the extra factor of dJ
619  Real fast_rhs_rho = -(temp_rhs_arr(i,j,k,0) + ( zflux_hi - zflux_lo ) * dzi);
620  Real fast_rhs_rhotheta = -( temp_rhs_arr(i,j,k,1) + myhalf *
621  ( zflux_hi * (prim(i,j,k) + prim(i,j,k+1))
622  - zflux_lo * (prim(i,j,k) + prim(i,j,k-1)) ) * dzi );
623 
624  cur_cons(i,j,k,0) *= (detJ_old(i,j,k)/detJ_new(i,j,k));
625  cur_cons(i,j,k,1) *= (detJ_old(i,j,k)/detJ_new(i,j,k));
626 
627  cur_cons(i,j,k,0) += dtau * ( slow_rhs_cons(i,j,k,0) + fast_rhs_rho / detJ_new(i,j,k));
628  cur_cons(i,j,k,1) += dtau * ( slow_rhs_cons(i,j,k,1) + fast_rhs_rhotheta / detJ_new(i,j,k));
629 
630  if (l_reflux) {
631  (flx_arr[2])(i,j,k,1) = (flx_arr[2])(i,j,k,0) * myhalf * (prim(i,j,k) + prim(i,j,k-1));
632  }
633 
634  if (k == vbx_hi.z) {
635  avg_zmom_arr(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
636  if (l_reflux) {
637  (flx_arr[2])(i,j,k+1,0) = zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
638  (flx_arr[2])(i,j,k+1,1) = (flx_arr[2])(i,j,k+1,0) * myhalf * (prim(i,j,k) + prim(i,j,k+1));
639  }
640  }
641 
642  // add in source terms for cell-centered conserved variables
643  cur_cons(i,j,k,Rho_comp) += dtau * cc_src_arr(i,j,k,Rho_comp);
644  cur_cons(i,j,k,RhoTheta_comp) += dtau * cc_src_arr(i,j,k,RhoTheta_comp);
645  });
646  } // end profile
647 
648  // We only add to the flux registers in the final RK step
649  if (l_reflux) {
650  int strt_comp_reflux = 0;
651  // For now we don't reflux (rho theta) because it seems to create issues at c/f boundaries
652  int num_comp_reflux = 1;
653  if (level < finest_level) {
654  fr_as_crse->CrseAdd(mfi,
655  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
656  dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
657  }
658  if (level > 0) {
659  fr_as_fine->FineAdd(mfi,
660  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
661  dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
662  }
663 
664  // This is necessary here so we don't go on to the next FArrayBox without
665  // having finished copying the fluxes into the FluxRegisters (since the fluxes
666  // are stored in temporary FArrayBox's)
667  Gpu::streamSynchronize();
668 
669  } // two-way coupling
670 
671  } // mfi
672  } // OMP
673 }
constexpr amrex::Real R_v
Definition: ERF_Constants.H:30
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
constexpr amrex::Real Gamma
Definition: ERF_Constants.H:38
@ 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
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:58
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:55
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
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
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real OmegaFromW(int &i, int &j, int &k, amrex::Real w, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:414
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:117
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:104
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:144
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:170
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int &i, int &j, int &k, amrex::Real omega, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:464
@ gpy
Definition: ERF_IndexDefines.H:185
@ gpx
Definition: ERF_IndexDefines.H:184
@ 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
@ qv
Definition: ERF_Kessler.H:29
@ q
Definition: ERF_WSM6.H:169
Here is the call graph for this function: