Function for computing the fast RHS with no terrain and variable vertical spacing
78 BL_PROFILE_REGION(
"erf_substep_S()");
80 const Box& domain = geom.Domain();
81 auto const domlo = lbound(domain);
82 auto const domhi = ubound(domain);
85 int ihi = domhi.x + 1;
87 int jhi = domhi.y + 1;
97 bool l_rayleigh_impl_for_w = (sinesq_stag_d !=
nullptr);
99 const Real*
dx = geom.CellSize();
100 const GpuArray<Real, AMREX_SPACEDIM>
dxInv = geom.InvCellSizeArray();
105 auto dz_ptr = stretched_dz_d.data();
108 const auto& dm = S_stage_data[
IntVars::cons].DistributionMap();
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));
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);
121 const Array<Real,AMREX_SPACEDIM> grav{
zero,
zero, -gravity};
122 const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
142 #pragma omp parallel if (Gpu::notInLaunchRegion())
144 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
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);
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);
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);
158 Box gbx = mfi.growntilebox(1);
159 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
164 theta_extrap(i,j,k) = prev_drho_theta(i,j,k);
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) );
172 theta_extrap(i,j,k) *= (
one + RvOverRd*
qv);
176 lagged_arr(i,j,k) = prev_drho_theta(i,j,k);
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);
193 #pragma omp parallel if (Gpu::notInLaunchRegion())
195 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
197 Box tbx = mfi.nodaltilebox(0);
198 Box tby = mfi.nodaltilebox(1);
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);
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);
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);
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);
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);
217 const Array4< Real>& avg_xmom_arr = avg_xmom.array(mfi);
218 const Array4< Real>& avg_ymom_arr = avg_ymom.array(mfi);
220 const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
222 const Array4<Real>& theta_extrap = extrap.array(mfi);
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);
231 if (nrk == 0 and step == 0) {
233 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
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;
239 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
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;
247 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
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;
254 Real q = (l_use_moisture) ?
myhalf * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) :
zero;
256 Real pi_c =
myhalf * (pi_stage_ca(i-1,j,k,0) + pi_stage_ca(i,j,k,0));
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);
263 avg_xmom_arr(i,j,k) += facinv*new_drho_u;
265 temp_cur_xmom_arr(i,j,k) = stage_xmom(i,j,k) + new_drho_u;
267 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
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;
274 Real q = (l_use_moisture) ?
myhalf * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) :
zero;
276 Real pi_c =
myhalf * (pi_stage_ca(i,j-1,k,0) + pi_stage_ca(i,j,k,0));
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);
283 avg_ymom_arr(i,j,k) += facinv*new_drho_v;
285 temp_cur_ymom_arr(i,j,k) = stage_ymom(i,j,k) + new_drho_v;
291 #pragma omp parallel if (Gpu::notInLaunchRegion())
294 std::array<FArrayBox,AMREX_SPACEDIM> flux;
297 Box bx = mfi.tilebox();
298 Box tbz = surroundingNodes(bx,2);
300 Box vbx = mfi.validbox();
301 const auto& vbx_hi = ubound(vbx);
303 const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
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);
310 const Array4<const Real>& prev_drho_theta = Delta_rho_theta.array(mfi);
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);
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);
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);
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);
325 const Array4< Real>& avg_zmom_arr = avg_zmom.array(mfi);
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);
334 RHS_fab.resize(tbz,1, The_Async_Arena());
337 soln_fab.resize(tbz,1, The_Async_Arena());
339 auto const& RHS_a = RHS_fab.array();
340 auto const& soln_a = soln_fab.array();
342 auto const& temp_rhs_arr = temp_rhs.array(mfi);
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);
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);
357 const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
358 flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
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);
367 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
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;
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));
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));
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));
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));
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);
409 ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
411 Real coeff_P = coeffP_a(i,j,k);
412 Real coeff_Q = coeffQ_a(i,j,k);
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);
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 );
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 );
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));
448 auto const lo = lbound(bx);
449 auto const hi = ubound(bx);
451 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
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);
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);
465 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
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);
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);
475 cur_zmom(i,j,hi.z+1) = stage_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
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);
483 for (
int j = lo.y; j <= hi.y; ++j) {
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);
489 for (
int k = lo.z+1; k <= hi.z+1; ++k) {
490 for (
int j = lo.y; j <= hi.y; ++j) {
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);
497 for (
int j = lo.y; j <= hi.y; ++j) {
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);
503 for (
int k = hi.z; k >= lo.z; --k) {
504 for (
int j = lo.y; j <= hi.y; ++j) {
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);
513 if (l_rayleigh_impl_for_w) {
514 ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
516 Real damping_coeff = l_damp_coef * dtau * sinesq_stag_d[k];
517 cur_zmom(i,j,k) /= (
one + damping_coeff);
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
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);
530 avg_zmom_arr(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
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));
537 avg_zmom_arr(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
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));
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)) );
552 int strt_comp_reflux = 0;
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);
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);
569 Gpu::streamSynchronize();
576 #pragma omp parallel if (Gpu::notInLaunchRegion())
578 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
580 const Box& bx = mfi.tilebox();
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);
589 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
601 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
615 const Array4<Real>& cur_xmom = S_data[
IntVars::xmom].array(mfi);
616 const Array4<Real>& cur_ymom = S_data[
IntVars::ymom].array(mfi);
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);
621 Box tbx = surroundingNodes(bx,0);
622 Box tby = surroundingNodes(bx,1);
625 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
627 cur_xmom(i,j,k) = temp_cur_xmom_arr(i,j,k);
629 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
631 cur_ymom(i,j,k) = temp_cur_ymom_arr(i,j,k);
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