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

Functions

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

Function Documentation

◆ erf_substep_NS()

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

Function for computing the fast RHS with no terrain and variable vertical spacing

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_previf step == 0, this is S_old, else the previous fast solution
[in]S_stage_datasolution at previous RK stage
[in]S_stage_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]stretched_dz_d
[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?
73 {
74  //
75  // NOTE: for step > 0, S_data and S_prev point to the same MultiFab data!!
76  //
77 
78  BL_PROFILE_REGION("erf_substep_S()");
79 
80  const Box& domain = geom.Domain();
81  auto const domlo = lbound(domain);
82  auto const domhi = ubound(domain);
83 
84  int ilo = domlo.x;
85  int ihi = domhi.x + 1;
86  int jlo = domlo.y;
87  int jhi = domhi.y + 1;
88 
89  Real beta_1 = myhalf * (one - beta_s); // multiplies explicit terms
90  Real beta_2 = myhalf * (one + beta_s); // multiplies implicit terms
91 
92  // How much do we project forward the (rho theta) that is used in the horizontal momentum equations
93  Real beta_d = Real(0.1);
94 
95  Real RvOverRd = R_v / R_d;
96 
97  bool l_rayleigh_impl_for_w = (sinesq_stag_d != nullptr);
98 
99  const Real* dx = geom.CellSize();
100  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
101 
102  Real dxi = dxInv[0];
103  Real dyi = dxInv[1];
104 
105  auto dz_ptr = stretched_dz_d.data();
106 
107  const auto& ba = S_stage_data[IntVars::cons].boxArray();
108  const auto& dm = S_stage_data[IntVars::cons].DistributionMap();
109 
110  MultiFab Delta_rho_theta( ba , dm, 1, 1);
111  MultiFab Delta_rho_w (convert(ba,IntVect(0,0,1)), dm, 1, IntVect(1,1,0));
112 
113  MultiFab coeff_A_mf(fast_coeffs, make_alias, 0, 1);
114  MultiFab inv_coeff_B_mf(fast_coeffs, make_alias, 1, 1);
115  MultiFab coeff_C_mf(fast_coeffs, make_alias, 2, 1);
116  MultiFab coeff_P_mf(fast_coeffs, make_alias, 3, 1);
117  MultiFab coeff_Q_mf(fast_coeffs, make_alias, 4, 1);
118 
119  // *************************************************************************
120  // Set gravity as a vector
121  const Array<Real,AMREX_SPACEDIM> grav{zero, zero, -gravity};
122  const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
123 
124  // This will hold theta extrapolated forward in time
125  MultiFab extrap(S_data[IntVars::cons].boxArray(),S_data[IntVars::cons].DistributionMap(),1,1);
126 
127  // This will hold the update for (rho) and (rho theta)
128  MultiFab temp_rhs(S_stage_data[IntVars::zmom].boxArray(),S_stage_data[IntVars::zmom].DistributionMap(),2,0);
129 
130  // This will hold the new x- and y-momenta temporarily (so that we don't overwrite values we need when tiling)
131  MultiFab temp_cur_xmom(S_stage_data[IntVars::xmom].boxArray(),S_stage_data[IntVars::xmom].DistributionMap(),1,0);
132  MultiFab temp_cur_ymom(S_stage_data[IntVars::ymom].boxArray(),S_stage_data[IntVars::ymom].DistributionMap(),1,0);
133 
134  // We assume that in the first step (nrk == 0) we are only doing one substep.
135  AMREX_ALWAYS_ASSERT(nrk > 0 || step == 0);
136 
137  // *************************************************************************
138  // First set up some arrays we'll need
139  // *************************************************************************
140 
141 #ifdef _OPENMP
142 #pragma omp parallel if (Gpu::notInLaunchRegion())
143 #endif
144  for ( MFIter mfi(S_stage_data[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
145  {
146  const Array4<const Real>& prev_cons = S_prev[IntVars::cons].const_array(mfi);
147  const Array4<const Real>& prev_zmom = S_prev[IntVars::zmom].const_array(mfi);
148 
149  const Array4<const Real>& stage_cons = S_stage_data[IntVars::cons].const_array(mfi);
150  const Array4<const Real>& stage_zmom = S_stage_data[IntVars::zmom].const_array(mfi);
151 
152  const Array4<Real>& prev_drho_w = Delta_rho_w.array(mfi);
153  const Array4<Real>& prev_drho_theta = Delta_rho_theta.array(mfi);
154  const Array4<Real>& lagged_arr = lagged_delta_rt.array(mfi);
155  const Array4<Real>& theta_extrap = extrap.array(mfi);
156  const Array4<const Real>& prim = S_stage_prim.const_array(mfi);
157 
158  Box gbx = mfi.growntilebox(1);
159  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
160  {
161  prev_drho_theta(i,j,k) = prev_cons(i,j,k,RhoTheta_comp) - stage_cons(i,j,k,RhoTheta_comp);
162 
163  if (step == 0) {
164  theta_extrap(i,j,k) = prev_drho_theta(i,j,k);
165  } else {
166  theta_extrap(i,j,k) = prev_drho_theta(i,j,k) + beta_d *
167  ( prev_drho_theta(i,j,k) - lagged_arr(i,j,k) );
168  }
169 
170  // NOTE: qv is not changing over the fast steps so we use the stage data
171  Real qv = (l_use_moisture) ? prim(i,j,k,PrimQ1_comp) : zero;
172  theta_extrap(i,j,k) *= (one + RvOverRd*qv);
173 
174  // We define lagged_delta_rt for our next step as the current delta_rt
175  // (after using it above to extrapolate theta for this step)
176  lagged_arr(i,j,k) = prev_drho_theta(i,j,k);
177  });
178 
179  // NOTE: We must do this here because for step > 0, prev_zmom and cur_zmom both point to the same data,
180  // so by the time we would use prev_zmom to define zflux, it would have already been over-written.
181  Box gtbz = mfi.nodaltilebox(2);
182  gtbz.grow(IntVect(1,1,0));
183  ParallelFor(gtbz, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
184  prev_drho_w(i,j,k) = prev_zmom(i,j,k) - stage_zmom(i,j,k);
185  });
186  } // mfi
187 
188  // *************************************************************************
189  // Define updates in the current RK stage
190  // *************************************************************************
191 
192 #ifdef _OPENMP
193 #pragma omp parallel if (Gpu::notInLaunchRegion())
194 #endif
195  for ( MFIter mfi(S_stage_data[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
196  {
197  Box tbx = mfi.nodaltilebox(0);
198  Box tby = mfi.nodaltilebox(1);
199 
200  const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
201  const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
202 
203  const Array4<const Real> & stage_xmom = S_stage_data[IntVars::xmom].const_array(mfi);
204  const Array4<const Real> & stage_ymom = S_stage_data[IntVars::ymom].const_array(mfi);
205  const Array4<const Real> & qt_arr = qt.const_array(mfi);
206 
207  const Array4<const Real>& slow_rhs_rho_u = S_slow_rhs[IntVars::xmom].const_array(mfi);
208  const Array4<const Real>& slow_rhs_rho_v = S_slow_rhs[IntVars::ymom].const_array(mfi);
209 
210  const Array4<Real>& temp_cur_xmom_arr = temp_cur_xmom.array(mfi);
211  const Array4<Real>& temp_cur_ymom_arr = temp_cur_ymom.array(mfi);
212 
213  const Array4<const Real>& prev_xmom = S_prev[IntVars::xmom].const_array(mfi);
214  const Array4<const Real>& prev_ymom = S_prev[IntVars::ymom].const_array(mfi);
215 
216  // These store the advection momenta which we will use to update the slow variables
217  const Array4< Real>& avg_xmom_arr = avg_xmom.array(mfi);
218  const Array4< Real>& avg_ymom_arr = avg_ymom.array(mfi);
219 
220  const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
221 
222  const Array4<Real>& theta_extrap = extrap.array(mfi);
223 
224  // Map factors
225  const Array4<const Real>& mf_ux = mapfac[MapFacType::u_x]->const_array(mfi);
226  const Array4<const Real>& mf_vy = mapfac[MapFacType::v_y]->const_array(mfi);
227 
228  // *********************************************************************
229  // Define updates in the RHS of {x, y, z}-momentum equations
230  // *********************************************************************
231  if (nrk == 0 and step == 0) { // prev == stage
232  ParallelFor(tbx, tby,
233  [=] AMREX_GPU_DEVICE (int i, int j, int k)
234  {
235  Real new_drho_u = dtau * slow_rhs_rho_u(i,j,k) + dtau * xmom_src_arr(i,j,k);;
236  avg_xmom_arr(i,j,k) += facinv*new_drho_u;
237  temp_cur_xmom_arr(i,j,k) = stage_xmom(i,j,k) + new_drho_u;
238  },
239  [=] AMREX_GPU_DEVICE (int i, int j, int k)
240  {
241  Real new_drho_v = dtau * slow_rhs_rho_v(i,j,k) + dtau * ymom_src_arr(i,j,k);
242  avg_ymom_arr(i,j,k) += facinv*new_drho_v;
243  temp_cur_ymom_arr(i,j,k) = stage_ymom(i,j,k) + new_drho_v;
244  });
245  } else {
246  ParallelFor(tbx, tby,
247  [=] AMREX_GPU_DEVICE (int i, int j, int k)
248  {
249  // Add (negative) gradient of (rho theta) multiplied by lagged "pi"
250  Real gpx = (l_real_bc && (level==0) && (i==ilo || i==ihi)) ? Real(0.) :
251  (theta_extrap(i,j,k) - theta_extrap(i-1,j,k))*dxi;
252  gpx *= mf_ux(i,j,0);
253 
254  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) : zero;
255 
256  Real pi_c = myhalf * (pi_stage_ca(i-1,j,k,0) + pi_stage_ca(i,j,k,0));
257  Real fast_rhs_rho_u = -Gamma * R_d * pi_c * gpx / (one + q);
258 
259  Real new_drho_u = prev_xmom(i,j,k) - stage_xmom(i,j,k)
260  + dtau * fast_rhs_rho_u + dtau * slow_rhs_rho_u(i,j,k)
261  + dtau * xmom_src_arr(i,j,k);
262 
263  avg_xmom_arr(i,j,k) += facinv*new_drho_u;
264 
265  temp_cur_xmom_arr(i,j,k) = stage_xmom(i,j,k) + new_drho_u;
266  },
267  [=] AMREX_GPU_DEVICE (int i, int j, int k)
268  {
269  // Add (negative) gradient of (rho theta) multiplied by lagged "pi"
270  Real gpy = (l_real_bc && (level==0) && (j==jlo || j==jhi)) ? Real(0.) :
271  (theta_extrap(i,j,k) - theta_extrap(i,j-1,k))*dyi;
272  gpy *= mf_vy(i,j,0);
273 
274  Real q = (l_use_moisture) ? myhalf * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) : zero;
275 
276  Real pi_c = myhalf * (pi_stage_ca(i,j-1,k,0) + pi_stage_ca(i,j,k,0));
277  Real fast_rhs_rho_v = -Gamma * R_d * pi_c * gpy / (one + q);
278 
279  Real new_drho_v = prev_ymom(i,j,k) - stage_ymom(i,j,k)
280  + dtau * fast_rhs_rho_v + dtau * slow_rhs_rho_v(i,j,k)
281  + dtau * ymom_src_arr(i,j,k);
282 
283  avg_ymom_arr(i,j,k) += facinv*new_drho_v;
284 
285  temp_cur_ymom_arr(i,j,k) = stage_ymom(i,j,k) + new_drho_v;
286  });
287  } // nrk > 0 and/or step > 0
288  } //mfi
289 
290 #ifdef _OPENMP
291 #pragma omp parallel if (Gpu::notInLaunchRegion())
292 #endif
293  {
294  std::array<FArrayBox,AMREX_SPACEDIM> flux;
295  for ( MFIter mfi(S_stage_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
296  {
297  Box bx = mfi.tilebox();
298  Box tbz = surroundingNodes(bx,2);
299 
300  Box vbx = mfi.validbox();
301  const auto& vbx_hi = ubound(vbx);
302 
303  const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
304 
305  const Array4<const Real>& stage_xmom = S_stage_data[IntVars::xmom].const_array(mfi);
306  const Array4<const Real>& stage_ymom = S_stage_data[IntVars::ymom].const_array(mfi);
307  const Array4<const Real>& stage_zmom = S_stage_data[IntVars::zmom].const_array(mfi);
308  const Array4<const Real> & prim = S_stage_prim.const_array(mfi);
309 
310  const Array4<const Real>& prev_drho_theta = Delta_rho_theta.array(mfi);
311 
312  const Array4<const Real>& prev_cons = S_prev[IntVars::cons].const_array(mfi);
313  const Array4<const Real>& stage_cons = S_stage_data[IntVars::cons].const_array(mfi);
314 
315  const Array4<const Real>& slow_rhs_cons = S_slow_rhs[IntVars::cons].const_array(mfi);
316  const Array4<const Real>& slow_rhs_rho_w = S_slow_rhs[IntVars::zmom].const_array(mfi);
317 
318  const Array4<const Real>& prev_zmom = S_prev[IntVars::zmom].const_array(mfi);
319  const Array4< Real>& cur_zmom = S_data[IntVars::zmom].array(mfi);
320 
321  const Array4<Real>& temp_cur_xmom_arr = temp_cur_xmom.array(mfi);
322  const Array4<Real>& temp_cur_ymom_arr = temp_cur_ymom.array(mfi);
323 
324  // These store the advection momenta which we will use to update the slow variables
325  const Array4< Real>& avg_zmom_arr = avg_zmom.array(mfi);
326 
327  // Map factors
328  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
329  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
330  const Array4<const Real>& mf_uy = mapfac[MapFacType::u_y]->const_array(mfi);
331  const Array4<const Real>& mf_vx = mapfac[MapFacType::v_x]->const_array(mfi);
332 
333  FArrayBox RHS_fab;
334  RHS_fab.resize(tbz,1, The_Async_Arena());
335 
336  FArrayBox soln_fab;
337  soln_fab.resize(tbz,1, The_Async_Arena());
338 
339  auto const& RHS_a = RHS_fab.array();
340  auto const& soln_a = soln_fab.array();
341 
342  auto const& temp_rhs_arr = temp_rhs.array(mfi);
343 
344  auto const& coeffA_a = coeff_A_mf.array(mfi);
345  auto const& inv_coeffB_a = inv_coeff_B_mf.array(mfi);
346  auto const& coeffC_a = coeff_C_mf.array(mfi);
347  auto const& coeffP_a = coeff_P_mf.array(mfi);
348  auto const& coeffQ_a = coeff_Q_mf.array(mfi);
349 
350  // *************************************************************************
351  // Define flux arrays for use in advection
352  // *************************************************************************
353  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
354  flux[dir].resize(surroundingNodes(bx,dir),2,The_Async_Arena());
355  flux[dir].setVal<RunOn::Device>(0);
356  }
357  const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
358  flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
359 
360  // *********************************************************************
361  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
362  Real xflux_lo = (temp_cur_xmom_arr(i ,j,k) - stage_xmom(i ,j,k)) / mf_uy(i ,j,0);
363  Real xflux_hi = (temp_cur_xmom_arr(i+1,j,k) - stage_xmom(i+1,j,k)) / mf_uy(i+1,j,0);
364  Real yflux_lo = (temp_cur_ymom_arr(i,j ,k) - stage_ymom(i,j ,k)) / mf_vx(i,j ,0);
365  Real yflux_hi = (temp_cur_ymom_arr(i,j+1,k) - stage_ymom(i,j+1,k)) / mf_vx(i,j+1,0);
366 
367  Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
368 
369  temp_rhs_arr(i,j,k,Rho_comp ) = ( xflux_hi - xflux_lo ) * dxi * mfsq
370  + ( yflux_hi - yflux_lo ) * dyi * mfsq;
371  temp_rhs_arr(i,j,k,RhoTheta_comp) = (( xflux_hi * (prim(i,j,k,0) + prim(i+1,j,k,0)) -
372  xflux_lo * (prim(i,j,k,0) + prim(i-1,j,k,0)) ) * dxi * mfsq +
373  ( yflux_hi * (prim(i,j,k,0) + prim(i,j+1,k,0)) -
374  yflux_lo * (prim(i,j,k,0) + prim(i,j-1,k,0)) ) * dyi * mfsq) * myhalf;
375 
376  if (l_reflux) {
377  (flx_arr[0])(i,j,k,0) = xflux_lo;
378  (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));
379 
380  (flx_arr[1])(i,j,k,0) = yflux_lo;
381  (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));
382 
383  if (i == vbx_hi.x) {
384  (flx_arr[0])(i+1,j,k,0) = xflux_hi;
385  (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));
386  }
387  if (j == vbx_hi.y) {
388  (flx_arr[1])(i,j+1,k,0) = yflux_hi;
389  (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));
390  }
391  }
392  });
393 
394  Box bx_shrunk_in_k = bx;
395  int klo = tbz.smallEnd(2);
396  int khi = tbz.bigEnd(2);
397  bx_shrunk_in_k.setSmall(2,klo+1);
398  bx_shrunk_in_k.setBig(2,khi-1);
399 
400  // Note that the notes use "g" to mean the magnitude of gravity, so it is positive
401  // We set grav_gpu[2] to be the vector component which is negative
402  // We define halfg to match the notes (which is why we take the absolute value)
403  Real halfg = std::abs(myhalf * grav_gpu[2]);
404 
405  // *********************************************************************
406  // fast_loop_on_shrunk
407  // *********************************************************************
408  //Note we don't act on the bottom or top boundaries of the domain
409  ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (int i, int j, int k)
410  {
411  Real coeff_P = coeffP_a(i,j,k);
412  Real coeff_Q = coeffQ_a(i,j,k);
413 
414  Real theta_t_lo = myhalf * ( prim(i,j,k-2,PrimTheta_comp) + prim(i,j,k-1,PrimTheta_comp) );
415  Real theta_t_mid = myhalf * ( prim(i,j,k-1,PrimTheta_comp) + prim(i,j,k ,PrimTheta_comp) );
416  Real theta_t_hi = myhalf * ( prim(i,j,k ,PrimTheta_comp) + prim(i,j,k+1,PrimTheta_comp) );
417 
418  Real Omega_kp1 = prev_zmom(i,j,k+1) - stage_zmom(i,j,k+1);
419  Real Omega_k = prev_zmom(i,j,k ) - stage_zmom(i,j,k );
420  Real Omega_km1 = prev_zmom(i,j,k-1) - stage_zmom(i,j,k-1);
421 
422  // line 2 last two terms (order dtau)
423  Real old_drho_k = prev_cons(i,j,k ,Rho_comp) - stage_cons(i,j,k ,Rho_comp);
424  Real old_drho_km1 = prev_cons(i,j,k-1,Rho_comp) - stage_cons(i,j,k-1,Rho_comp);
425  Real R0_tmp = coeff_P * prev_drho_theta(i,j,k) + coeff_Q * prev_drho_theta(i,j,k-1)
426  - halfg * ( old_drho_k + old_drho_km1 );
427 
428  // lines 3-5 residuals (order dtau^2) one <-> beta_2
429  Real R1_tmp = halfg * (-slow_rhs_cons(i,j,k ,Rho_comp) - slow_rhs_cons(i,j,k-1,Rho_comp)
430  + temp_rhs_arr(i,j,k ,Rho_comp) + temp_rhs_arr(i,j,k-1,Rho_comp) )
431  + ( coeff_P * (slow_rhs_cons(i,j,k ,RhoTheta_comp) - temp_rhs_arr(i,j,k ,RhoTheta_comp)) +
432  coeff_Q * (slow_rhs_cons(i,j,k-1,RhoTheta_comp) - temp_rhs_arr(i,j,k-1,RhoTheta_comp)) );
433 
434  // lines 6&7 consolidated (reuse Omega & metrics) (order dtau^2)
435  Real dz_inv = one / dz_ptr[k];
436  R1_tmp += beta_1 * dz_inv * ( (Omega_kp1 - Omega_km1) * halfg
437  -(Omega_kp1*theta_t_hi - Omega_k *theta_t_mid) * coeff_P
438  -(Omega_k *theta_t_mid - Omega_km1*theta_t_lo ) * coeff_Q );
439 
440  // line 1
441  RHS_a(i,j,k) = Omega_k + dtau * (slow_rhs_rho_w(i,j,k) + R0_tmp + dtau * beta_2 * R1_tmp + zmom_src_arr(i,j,k));
442 
443  }); // bx_shrunk_in_k
444 
445  Box b2d = tbz; // Copy constructor
446  b2d.setRange(2,0);
447 
448  auto const lo = lbound(bx);
449  auto const hi = ubound(bx);
450 
451  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int)
452  {
453  // w at bottom boundary of grid is 0 if at domain boundary, otherwise w = w_old + dtau * slow_rhs
454  RHS_a (i,j,lo.z) = prev_zmom(i,j,lo.z) - stage_zmom(i,j,lo.z)
455  + dtau * slow_rhs_rho_w(i,j,lo.z)
456  + dtau * zmom_src_arr(i,j,lo.z);
457 
458  // w at top boundary of grid is 0 if at domain boundary, otherwise w = w_old + dtau * slow_rhs
459  RHS_a (i,j,hi.z+1) = prev_zmom(i,j,hi.z+1) - stage_zmom(i,j,hi.z+1)
460  + dtau * slow_rhs_rho_w(i,j,hi.z+1)
461  + dtau * zmom_src_arr(i,j,hi.z+1);
462  }); // b2d
463 
464 #ifdef AMREX_USE_GPU
465  ParallelFor(b2d, [=] AMREX_GPU_DEVICE (int i, int j, int)
466  {
467  // w = specified Dirichlet value at k = lo.z
468  soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
469  cur_zmom(i,j,lo.z) = stage_zmom(i,j,lo.z) + soln_a(i,j,lo.z);
470 
471  for (int k = lo.z+1; k <= hi.z+1; k++) {
472  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);
473  }
474 
475  cur_zmom(i,j,hi.z+1) = stage_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
476 
477  for (int k = hi.z; k >= lo.z; k--) {
478  soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) *soln_a(i,j,k+1);
479  cur_zmom(i,j,k) = stage_zmom(i,j,k) + soln_a(i,j,k);
480  }
481  }); // b2d
482 #else
483  for (int j = lo.y; j <= hi.y; ++j) {
484  AMREX_PRAGMA_SIMD
485  for (int i = lo.x; i <= hi.x; ++i) {
486  soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
487  }
488  }
489  for (int k = lo.z+1; k <= hi.z+1; ++k) {
490  for (int j = lo.y; j <= hi.y; ++j) {
491  AMREX_PRAGMA_SIMD
492  for (int i = lo.x; i <= hi.x; ++i) {
493  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);
494  }
495  }
496  }
497  for (int j = lo.y; j <= hi.y; ++j) {
498  AMREX_PRAGMA_SIMD
499  for (int i = lo.x; i <= hi.x; ++i) {
500  cur_zmom(i,j,hi.z+1) = stage_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
501  }
502  }
503  for (int k = hi.z; k >= lo.z; --k) {
504  for (int j = lo.y; j <= hi.y; ++j) {
505  AMREX_PRAGMA_SIMD
506  for (int i = lo.x; i <= hi.x; ++i) {
507  soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
508  cur_zmom(i,j,k) = stage_zmom(i,j,k) + soln_a(i,j,k);
509  }
510  }
511  }
512 #endif
513  if (l_rayleigh_impl_for_w) {
514  ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (int i, int j, int k)
515  {
516  Real damping_coeff = l_damp_coef * dtau * sinesq_stag_d[k];
517  cur_zmom(i,j,k) /= (one + damping_coeff);
518  });
519  }
520 
521  // **************************************************************************
522  // Define updates in the RHS of rho and (rho theta)
523  // **************************************************************************
524  const Array4<Real>& prev_drho_w = Delta_rho_w.array(mfi);
525  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
526  {
527  Real zflux_lo = beta_2 * soln_a(i,j,k ) + beta_1 * prev_drho_w(i,j,k );
528  Real zflux_hi = beta_2 * soln_a(i,j,k+1) + beta_1 * prev_drho_w(i,j,k+1);
529 
530  avg_zmom_arr(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
531  if (l_reflux) {
532  (flx_arr[2])(i,j,k,0) = zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
533  (flx_arr[2])(i,j,k,1) = (flx_arr[2])(i,j,k,0) * myhalf * (prim(i,j,k) + prim(i,j,k-1));
534  }
535 
536  if (k == vbx_hi.z) {
537  avg_zmom_arr(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
538  if (l_reflux) {
539  (flx_arr[2])(i,j,k+1,0) = zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
540  (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));
541  }
542  }
543 
544  Real dz_inv = one / dz_ptr[k];
545  temp_rhs_arr(i,j,k,Rho_comp ) += dz_inv * ( zflux_hi - zflux_lo );
546  temp_rhs_arr(i,j,k,RhoTheta_comp) += myhalf * dz_inv * ( zflux_hi * (prim(i,j,k) + prim(i,j,k+1))
547  - zflux_lo * (prim(i,j,k) + prim(i,j,k-1)) );
548  });
549 
550  // We only add to the flux registers in the final RK step
551  if (l_reflux) {
552  int strt_comp_reflux = 0;
553  // For now we don't reflux (rho theta) because it seems to create issues at c/f boundaries
554  int num_comp_reflux = 1;
555  if (level < finest_level) {
556  fr_as_crse->CrseAdd(mfi,
557  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
558  dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
559  }
560  if (level > 0) {
561  fr_as_fine->FineAdd(mfi,
562  {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
563  dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
564  }
565 
566  // This is necessary here so we don't go on to the next FArrayBox without
567  // having finished copying the fluxes into the FluxRegisters (since the fluxes
568  // are stored in temporary FArrayBox's)
569  Gpu::streamSynchronize();
570 
571  } // two-way coupling
572  } // mfi
573  } // OMP
574 
575 #ifdef _OPENMP
576 #pragma omp parallel if (Gpu::notInLaunchRegion())
577 #endif
578  for ( MFIter mfi(S_stage_data[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
579  {
580  const Box& bx = mfi.tilebox();
581 
582  const Array4< Real>& cur_cons = S_data[IntVars::cons].array(mfi);
583  const Array4<const Real>& prev_cons = S_prev[IntVars::cons].const_array(mfi);
584  auto const& temp_rhs_arr = temp_rhs.const_array(mfi);
585  auto const& slow_rhs_cons = S_slow_rhs[IntVars::cons].const_array(mfi);
586  const Array4<Real const>& cc_src_arr = cc_src.const_array(mfi);
587 
588  if (step == 0) {
589  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
590  {
591  cur_cons(i,j,k,Rho_comp) = prev_cons(i,j,k,Rho_comp) +
592  dtau * (slow_rhs_cons(i,j,k,Rho_comp) - temp_rhs_arr(i,j,k,Rho_comp));
593  cur_cons(i,j,k,RhoTheta_comp) = prev_cons(i,j,k,RhoTheta_comp) +
594  dtau * (slow_rhs_cons(i,j,k,RhoTheta_comp) - temp_rhs_arr(i,j,k,RhoTheta_comp));
595 
596  // add in source terms for cell-centered conserved variables
597  cur_cons(i,j,k,Rho_comp) += dtau * cc_src_arr(i,j,k,Rho_comp);
598  cur_cons(i,j,k,RhoTheta_comp) += dtau * cc_src_arr(i,j,k,RhoTheta_comp);
599  });
600  } else {
601  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
602  {
603  //
604  // We didn't need to set cur_cons = prev_cons above because they point to the same data for step > 0
605  //
606  cur_cons(i,j,k,Rho_comp) += dtau * (slow_rhs_cons(i,j,k,Rho_comp) - temp_rhs_arr(i,j,k,Rho_comp));
607  cur_cons(i,j,k,RhoTheta_comp) += dtau * (slow_rhs_cons(i,j,k,RhoTheta_comp) - temp_rhs_arr(i,j,k,RhoTheta_comp));
608 
609  // add in source terms for cell-centered conserved variables
610  cur_cons(i,j,k,Rho_comp) += dtau * cc_src_arr(i,j,k,Rho_comp);
611  cur_cons(i,j,k,RhoTheta_comp) += dtau * cc_src_arr(i,j,k,RhoTheta_comp);
612  });
613  } // step = 0
614 
615  const Array4<Real>& cur_xmom = S_data[IntVars::xmom].array(mfi);
616  const Array4<Real>& cur_ymom = S_data[IntVars::ymom].array(mfi);
617 
618  const Array4<Real const>& temp_cur_xmom_arr = temp_cur_xmom.const_array(mfi);
619  const Array4<Real const>& temp_cur_ymom_arr = temp_cur_ymom.const_array(mfi);
620 
621  Box tbx = surroundingNodes(bx,0);
622  Box tby = surroundingNodes(bx,1);
623 
624  ParallelFor(tbx, tby,
625  [=] AMREX_GPU_DEVICE (int i, int j, int k)
626  {
627  cur_xmom(i,j,k) = temp_cur_xmom_arr(i,j,k);
628  },
629  [=] AMREX_GPU_DEVICE (int i, int j, int k)
630  {
631  cur_ymom(i,j,k) = temp_cur_ymom_arr(i,j,k);
632  });
633 
634  } // mfi
635 }
constexpr amrex::Real R_v
Definition: ERF_Constants.H:30
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
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
constexpr amrex::Real Gamma
Definition: ERF_Constants.H:38
@ v_x
Definition: ERF_DataStruct.H:24
@ u_y
Definition: ERF_DataStruct.H:25
@ v_y
Definition: ERF_DataStruct.H:25
@ m_y
Definition: ERF_DataStruct.H:25
@ u_x
Definition: ERF_DataStruct.H:24
@ m_x
Definition: ERF_DataStruct.H:24
#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
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
@ 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: