81 BL_PROFILE_REGION(
"erf_fast_rhs_MT()");
83 Real beta_1 = 0.5 * (1.0 - beta_s);
84 Real beta_2 = 0.5 * (1.0 + beta_s);
89 const Real* dx = geom.CellSize();
90 const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
96 MultiFab coeff_A_mf(fast_coeffs, make_alias, 0, 1);
97 MultiFab inv_coeff_B_mf(fast_coeffs, make_alias, 1, 1);
98 MultiFab coeff_C_mf(fast_coeffs, make_alias, 2, 1);
99 MultiFab coeff_P_mf(fast_coeffs, make_alias, 3, 1);
100 MultiFab coeff_Q_mf(fast_coeffs, make_alias, 4, 1);
104 const Array<Real,AMREX_SPACEDIM> grav{0.0, 0.0, -gravity};
105 const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
115 #pragma omp parallel if (Gpu::notInLaunchRegion())
118 FArrayBox temp_rhs_fab;
123 std::array<FArrayBox,AMREX_SPACEDIM> flux;
127 for ( MFIter mfi(S_stg_data[
IntVars::cons],
false); mfi.isValid(); ++mfi)
129 Box bx = mfi.tilebox();
130 Box tbx = surroundingNodes(bx,0);
131 Box tby = surroundingNodes(bx,1);
132 Box tbz = surroundingNodes(bx,2);
134 Box vbx = mfi.validbox();
135 const auto& vbx_hi = ubound(vbx);
137 const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
138 const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
139 const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
140 const Array4<Real const>& cc_src_arr = cc_src.const_array(mfi);
142 const Array4<const Real> & stg_cons = S_stg_data[
IntVars::cons].const_array(mfi);
143 const Array4<const Real> & stg_xmom = S_stg_data[
IntVars::xmom].const_array(mfi);
144 const Array4<const Real> & stg_ymom = S_stg_data[
IntVars::ymom].const_array(mfi);
145 const Array4<const Real> & stg_zmom = S_stg_data[
IntVars::zmom].const_array(mfi);
146 const Array4<const Real> & prim = S_stg_prim.const_array(mfi);
148 const Array4<const Real>& slow_rhs_cons = S_slow_rhs[
IntVars::cons].const_array(mfi);
149 const Array4<const Real>& slow_rhs_rho_u = S_slow_rhs[
IntVars::xmom].const_array(mfi);
150 const Array4<const Real>& slow_rhs_rho_v = S_slow_rhs[
IntVars::ymom].const_array(mfi);
151 const Array4<const Real>& slow_rhs_rho_w = S_slow_rhs[
IntVars::zmom].const_array(mfi);
153 const Array4<Real>& cur_cons = S_data[
IntVars::cons].array(mfi);
154 const Array4<Real>& cur_xmom = S_data[
IntVars::xmom].array(mfi);
155 const Array4<Real>& cur_ymom = S_data[
IntVars::ymom].array(mfi);
156 const Array4<Real>& cur_zmom = S_data[
IntVars::zmom].array(mfi);
158 const Array4<Real>& lagged_delta_rt = S_scratch[
IntVars::cons].array(mfi);
160 const Array4<const Real>& prev_cons = S_prev[
IntVars::cons].const_array(mfi);
161 const Array4<const Real>& prev_xmom = S_prev[
IntVars::xmom].const_array(mfi);
162 const Array4<const Real>& prev_ymom = S_prev[
IntVars::ymom].const_array(mfi);
163 const Array4<const Real>& prev_zmom = S_prev[
IntVars::zmom].const_array(mfi);
166 const Array4<Real>& avg_xmom = S_scratch[
IntVars::xmom].array(mfi);
167 const Array4<Real>& avg_ymom = S_scratch[
IntVars::ymom].array(mfi);
168 const Array4<Real>& avg_zmom = S_scratch[
IntVars::zmom].array(mfi);
170 const Array4<const Real>& z_nd_old = z_phys_nd_old->const_array(mfi);
171 const Array4<const Real>& z_nd_new = z_phys_nd_new->const_array(mfi);
172 const Array4<const Real>& z_nd_stg = z_phys_nd_stg->const_array(mfi);
173 const Array4<const Real>& detJ_old = detJ_cc_old->const_array(mfi);
174 const Array4<const Real>& detJ_new = detJ_cc_new->const_array(mfi);
175 const Array4<const Real>& detJ_stg = detJ_cc_stg->const_array(mfi);
177 const Array4<const Real>& z_t_arr = z_t_rk->const_array(mfi);
178 const Array4<const Real>& zp_t_arr = z_t_pert->const_array(mfi);
180 const Array4< Real>& omega_arr = Omega.array(mfi);
182 const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
184 const Array4<Real>& theta_extrap = extrap.array(mfi);
187 const Array4<const Real>& mf_mx = mapfac[
MapFacType::m_x]->const_array(mfi);
188 const Array4<const Real>& mf_my = mapfac[
MapFacType::m_y]->const_array(mfi);
189 const Array4<const Real>& mf_ux = mapfac[
MapFacType::u_x]->const_array(mfi);
190 const Array4<const Real>& mf_vy = mapfac[
MapFacType::v_y]->const_array(mfi);
194 Box gbx = mfi.tilebox(); gbx.grow(1);
195 Box gtbx = mfi.nodaltilebox(0); gtbx.grow(1); gtbx.setSmall(2,0);
196 Box gtby = mfi.nodaltilebox(1); gtby.grow(1); gtby.setSmall(2,0);
199 BL_PROFILE(
"fast_rhs_copies_0");
202 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
207 theta_extrap(i,j,k) = delta_rt;
212 }
else if (use_lagged_delta_rt) {
215 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
217 theta_extrap(i,j,k) = delta_rt + beta_d * (delta_rt - lagged_delta_rt(i,j,k,
RhoTheta_comp));
225 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
231 RHS_fab.resize (tbz,1, The_Async_Arena());
232 soln_fab.resize (tbz,1, The_Async_Arena());
233 temp_rhs_fab.resize(tbz,2, The_Async_Arena());
235 auto const& RHS_a = RHS_fab.array();
236 auto const& soln_a = soln_fab.array();
237 auto const& temp_rhs_arr = temp_rhs_fab.array();
239 auto const& coeffA_a = coeff_A_mf.array(mfi);
240 auto const& inv_coeffB_a = inv_coeff_B_mf.array(mfi);
241 auto const& coeffC_a = coeff_C_mf.array(mfi);
242 auto const& coeffP_a = coeff_P_mf.array(mfi);
243 auto const& coeffQ_a = coeff_Q_mf.array(mfi);
249 BL_PROFILE(
"fast_rhs_xymom_T");
250 ParallelFor(tbx, tby,
251 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
256 Real gp_xi = (theta_extrap(i,j,k) - theta_extrap(i-1,j,k)) * dxi;
257 Real gp_zeta_on_iface = (k == 0) ?
258 0.5 * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
259 -theta_extrap(i-1,j,k ) - theta_extrap(i,j,k ) ) :
260 0.25 * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
261 -theta_extrap(i-1,j,k-1) - theta_extrap(i,j,k-1) );
262 Real
gpx = h_zeta_old * gp_xi - h_xi_old * gp_zeta_on_iface;
265 if (l_use_moisture) {
271 Real pi_c = 0.5 * (pi_stage_ca(i-1,j,k,0) + pi_stage_ca(i ,j,k,0));
275 cur_xmom(i,j,k) = h_zeta_old * prev_xmom(i,j,k) + dtau * fast_rhs_rho_u
276 + dtau * slow_rhs_rho_u(i,j,k)
277 + dtau * xmom_src_arr(i,j,k);
279 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
284 Real gp_eta = (theta_extrap(i,j,k) -theta_extrap(i,j-1,k)) * dyi;
285 Real gp_zeta_on_jface = (k == 0) ?
286 0.5 * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
287 -theta_extrap(i,j,k ) - theta_extrap(i,j-1,k ) ) :
288 0.25 * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
289 -theta_extrap(i,j,k-1) - theta_extrap(i,j-1,k-1) );
290 Real
gpy = h_zeta_old * gp_eta - h_eta_old * gp_zeta_on_jface;
293 if (l_use_moisture) {
299 Real pi_c = 0.5 * (pi_stage_ca(i,j-1,k,0) + pi_stage_ca(i,j ,k,0));
303 cur_ymom(i, j, k) = h_zeta_old * prev_ymom(i,j,k) + dtau * fast_rhs_rho_v
304 + dtau * slow_rhs_rho_v(i,j,k)
305 + dtau * ymom_src_arr(i,j,k);
312 for (
int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
313 flux[dir].resize(surroundingNodes(bx,dir),2);
314 flux[dir].setVal<RunOn::Device>(0.);
316 const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
317 flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
321 BL_PROFILE(
"fast_T_making_rho_rhs");
322 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
326 Real xflux_lo = cur_xmom(i ,j,k) - stg_xmom(i ,j,k)*h_zeta_stg_xlo;
327 Real xflux_hi = cur_xmom(i+1,j,k) - stg_xmom(i+1,j,k)*h_zeta_stg_xhi;
331 Real yflux_lo = cur_ymom(i,j ,k) - stg_ymom(i,j ,k)*h_zeta_stg_ylo;
332 Real yflux_hi = cur_ymom(i,j+1,k) - stg_ymom(i,j+1,k)*h_zeta_stg_yhi;
335 temp_rhs_arr(i,j,k,0) = ( xflux_hi - xflux_lo ) * dxi + ( yflux_hi - yflux_lo ) * dyi;
336 temp_rhs_arr(i,j,k,1) = (( xflux_hi * (prim(i,j,k,0) + prim(i+1,j,k,0)) -
337 xflux_lo * (prim(i,j,k,0) + prim(i-1,j,k,0)) ) * dxi +
338 ( yflux_hi * (prim(i,j,k,0) + prim(i,j+1,k,0)) -
339 yflux_lo * (prim(i,j,k,0) + prim(i,j-1,k,0)) ) * dyi) * 0.5;
342 (flx_arr[0])(i,j,k,0) = xflux_lo;
343 (flx_arr[0])(i,j,k,1) = (flx_arr[0])(i ,j,k,0) * 0.5 * (prim(i,j,k,0) + prim(i-1,j,k,0));
345 (flx_arr[1])(i,j,k,0) = yflux_lo;
346 (flx_arr[1])(i,j,k,1) = (flx_arr[1])(i,j ,k,0) * 0.5 * (prim(i,j,k,0) + prim(i,j-1,k,0));
349 (flx_arr[0])(i+1,j,k,0) = xflux_hi;
350 (flx_arr[0])(i+1,j,k,1) = (flx_arr[0])(i+1,j,k,0) * 0.5 * (prim(i,j,k,0) + prim(i+1,j,k,0));
353 (flx_arr[1])(i,j+1,k,0) = yflux_hi;
354 (flx_arr[1])(i,j+1,k,1) = (flx_arr[1])(i,j+1,k,0) * 0.5 * (prim(i,j,k,0) + prim(i,j+1,k,0));
364 Box gbxo = mfi.nodaltilebox(2);
366 BL_PROFILE(
"fast_MT_making_omega");
367 Box gbxo_lo = gbxo; gbxo_lo.setBig(2,0);
368 ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
369 omega_arr(i,j,k) = 0.;
371 Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
372 ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
373 omega_arr(i,j,k) = prev_zmom(i,j,k) - stg_zmom(i,j,k) - zp_t_arr(i,j,k);
375 Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
376 ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
378 (
OmegaFromW(i,j,k,prev_zmom(i,j,k),prev_xmom,prev_ymom,mf_ux,mf_vy,z_nd_old,dxInv)
379 -
OmegaFromW(i,j,k, stg_zmom(i,j,k), stg_xmom, stg_ymom,mf_ux,mf_vy,z_nd_old,dxInv) )
385 ParallelFor(tbx, tby,
386 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
389 cur_xmom(i, j, k) /= h_zeta_new;
390 avg_xmom(i,j,k) += facinv*(cur_xmom(i,j,k) - stg_xmom(i,j,k));
392 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
395 cur_ymom(i, j, k) /= h_zeta_new;
396 avg_ymom(i,j,k) += facinv*(cur_ymom(i,j,k) - stg_ymom(i,j,k));
399 Box bx_shrunk_in_k = bx;
400 int klo = tbz.smallEnd(2);
401 int khi = tbz.bigEnd(2);
402 bx_shrunk_in_k.setSmall(2,klo+1);
403 bx_shrunk_in_k.setBig(2,khi-1);
408 Real halfg = std::abs(0.5 * grav_gpu[2]);
411 BL_PROFILE(
"fast_loop_on_shrunk_t");
413 ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
415 Real dJ_old_kface = 0.5 * (detJ_old(i,j,k) + detJ_old(i,j,k-1));
416 Real dJ_new_kface = 0.5 * (detJ_new(i,j,k) + detJ_new(i,j,k-1));
417 Real dJ_stg_kface = 0.5 * (detJ_stg(i,j,k) + detJ_stg(i,j,k-1));
419 Real coeff_P = coeffP_a(i,j,k);
420 Real coeff_Q = coeffQ_a(i,j,k);
422 if (l_use_moisture) {
425 coeff_P /= (1.0 + q);
426 coeff_Q /= (1.0 + q);
434 Real R0_tmp = coeff_P * cur_cons(i,j,k ,
RhoTheta_comp) * dJ_old_kface
438 - halfg * ( cur_cons(i,j,k,
Rho_comp) + cur_cons(i,j,k-1,
Rho_comp) ) * dJ_old_kface
439 + halfg * ( stg_cons(i,j,k,
Rho_comp) + stg_cons(i,j,k-1,
Rho_comp) ) * dJ_stg_kface;
442 Real R1_tmp = - halfg * ( slow_rhs_cons(i,j,k ,
Rho_comp) +
447 Real Omega_kp1 = omega_arr(i,j,k+1);
448 Real Omega_k = omega_arr(i,j,k );
449 Real Omega_km1 = omega_arr(i,j,k-1);
452 R1_tmp += ( halfg ) *
453 ( beta_1 * dzi * (Omega_kp1 - Omega_km1) + temp_rhs_arr(i,j,k,
Rho_comp) + temp_rhs_arr(i,j,k-1,
Rho_comp));
457 coeff_P * ( beta_1 * dzi * (Omega_kp1*theta_t_hi - Omega_k*theta_t_mid) + temp_rhs_arr(i,j,k ,
RhoTheta_comp) ) +
458 coeff_Q * ( beta_1 * dzi * (Omega_k*theta_t_mid - Omega_km1*theta_t_lo) + temp_rhs_arr(i,j,k-1,
RhoTheta_comp) ) );
461 RHS_a(i,j,k) = dJ_old_kface * prev_zmom(i,j,k) - dJ_stg_kface * stg_zmom(i,j,k)
462 + dtau * (slow_rhs_rho_w(i,j,k) + R0_tmp + dtau*beta_2*R1_tmp + zmom_src_arr(i,j,k));
465 Real UppVpp = dJ_new_kface *
OmegaFromW(i,j,k,0.,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
466 -dJ_stg_kface *
OmegaFromW(i,j,k,0.,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
467 RHS_a(i,j,k) += UppVpp;
474 auto const lo = lbound(bx);
475 auto const hi = ubound(bx);
478 BL_PROFILE(
"fast_rhs_b2d_loop_t");
481 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
484 Real rho_on_bdy = 0.5 * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
485 RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,0);
487 soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
490 RHS_a(i,j,hi.z+1) = 0.0;
492 for (
int k = lo.z+1; k <= hi.z+1; k++) {
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);
496 for (
int k = hi.z; k >= lo.z; k--) {
497 soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
501 cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
504 for (
int j = lo.y; j <= hi.y; ++j) {
506 for (
int i = lo.x; i <= hi.x; ++i) {
508 Real rho_on_bdy = 0.5 * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
509 RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,lo.z);
511 soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
515 for (
int j = lo.y; j <= hi.y; ++j) {
517 for (
int i = lo.x; i <= hi.x; ++i) {
518 RHS_a (i,j,hi.z+1) = 0.0;
521 for (
int k = lo.z+1; k <= hi.z+1; ++k) {
522 for (
int j = lo.y; j <= hi.y; ++j) {
524 for (
int i = lo.x; i <= hi.x; ++i) {
525 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);
529 for (
int k = hi.z; k >= lo.z; --k) {
530 for (
int j = lo.y; j <= hi.y; ++j) {
532 for (
int i = lo.x; i <= hi.x; ++i) {
533 soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
539 for (
int j = lo.y; j <= hi.y; ++j) {
541 for (
int i = lo.x; i <= hi.x; ++i) {
542 cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
549 BL_PROFILE(
"fast_rhs_new_drhow");
551 ParallelFor(tbz, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
553 Real rho_on_face = 0.5 * (cur_cons(i,j,k,
Rho_comp) + cur_cons(i,j,k-1,
Rho_comp));
556 cur_zmom(i,j,k) =
WFromOmega(i,j,k,rho_on_face*(z_t_arr(i,j,k)+zp_t_arr(i,j,k)),
557 cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv);
564 Real UppVpp =
WFromOmega(i,j,k,0.0,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
565 -
WFromOmega(i,j,k,0.0,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
566 Real wpp = soln_a(i,j,k) + UppVpp;
567 Real dJ_old_kface = 0.5 * (detJ_old(i,j,k) + detJ_old(i,j,k-1));
568 Real dJ_new_kface = 0.5 * (detJ_new(i,j,k) + detJ_new(i,j,k-1));
570 cur_zmom(i,j,k) = dJ_old_kface * (stg_zmom(i,j,k) + wpp);
571 cur_zmom(i,j,k) /= dJ_new_kface;
573 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)
574 -
OmegaFromW(i,j,k,stg_zmom(i,j,k),stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
575 soln_a(i,j,k) -= rho_on_face * zp_t_arr(i,j,k);
584 BL_PROFILE(
"fast_rho_final_update");
585 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
587 Real zflux_lo = beta_2 * soln_a(i,j,k ) + beta_1 * omega_arr(i,j,k);
588 Real zflux_hi = beta_2 * soln_a(i,j,k+1) + beta_1 * omega_arr(i,j,k+1);
592 avg_zmom(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
594 (flx_arr[2])(i,j,k,0) = zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
600 Real fast_rhs_rho = -(temp_rhs_arr(i,j,k,0) + ( zflux_hi - zflux_lo ) * dzi);
601 Real temp_rho = detJ_old(i,j,k) * cur_cons(i,j,k,0) +
602 dtau * ( slow_rhs_cons(i,j,k,0) + fast_rhs_rho );
603 cur_cons(i,j,k,0) = temp_rho / detJ_new(i,j,k);
608 Real fast_rhs_rhotheta = -( temp_rhs_arr(i,j,k,1) + 0.5 *
609 ( zflux_hi * (prim(i,j,k) + prim(i,j,k+1))
610 - zflux_lo * (prim(i,j,k) + prim(i,j,k-1)) ) * dzi );
611 Real temp_rth = detJ_old(i,j,k) * cur_cons(i,j,k,1) +
612 dtau * ( slow_rhs_cons(i,j,k,1) + fast_rhs_rhotheta );
613 cur_cons(i,j,k,1) = temp_rth / detJ_new(i,j,k);
615 (flx_arr[2])(i,j,k,1) = (flx_arr[2])(i,j,k,0) * 0.5 * (prim(i,j,k) + prim(i,j,k-1));
619 avg_zmom(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
621 (flx_arr[2])(i,j,k+1,0) = zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
622 (flx_arr[2])(i,j,k+1,1) = (flx_arr[2])(i,j,k+1,0) * 0.5 * (prim(i,j,k) + prim(i,j,k+1));
634 int strt_comp_reflux = 0;
635 int num_comp_reflux = 2;
636 if (level < finest_level) {
637 fr_as_crse->CrseAdd(mfi,
638 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
639 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
642 fr_as_fine->FineAdd(mfi,
643 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
644 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
650 Gpu::streamSynchronize();
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
constexpr amrex::Real Gamma
Definition: ERF_Constants.H:19
@ v_y
Definition: ERF_DataStruct.H:22
@ m_y
Definition: ERF_DataStruct.H:22
@ u_x
Definition: ERF_DataStruct.H:21
@ m_x
Definition: ERF_DataStruct.H:21
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:53
#define PrimQ2_comp
Definition: ERF_IndexDefines.H:54
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
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:415
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:110
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:96
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:139
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:167
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:465
@ gpy
Definition: ERF_IndexDefines.H:151
@ gpx
Definition: ERF_IndexDefines.H:150
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159