86 BL_PROFILE_REGION(
"erf_substep_MT()");
88 Real beta_1 = 0.5 * (1.0 - beta_s);
89 Real beta_2 = 0.5 * (1.0 + beta_s);
96 const Real* dx = geom.CellSize();
97 const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
103 MultiFab coeff_A_mf(fast_coeffs, make_alias, 0, 1);
104 MultiFab inv_coeff_B_mf(fast_coeffs, make_alias, 1, 1);
105 MultiFab coeff_C_mf(fast_coeffs, make_alias, 2, 1);
106 MultiFab coeff_P_mf(fast_coeffs, make_alias, 3, 1);
107 MultiFab coeff_Q_mf(fast_coeffs, make_alias, 4, 1);
111 const Array<Real,AMREX_SPACEDIM> grav{0.0, 0.0, -gravity};
112 const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
122 #pragma omp parallel if (Gpu::notInLaunchRegion())
125 FArrayBox temp_rhs_fab;
130 std::array<FArrayBox,AMREX_SPACEDIM> flux;
134 for ( MFIter mfi(S_stg_data[
IntVars::cons],
false); mfi.isValid(); ++mfi)
136 Box bx = mfi.tilebox();
137 Box tbx = surroundingNodes(bx,0);
138 Box tby = surroundingNodes(bx,1);
139 Box tbz = surroundingNodes(bx,2);
141 Box vbx = mfi.validbox();
142 const auto& vbx_hi = ubound(vbx);
144 const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
145 const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
146 const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
147 const Array4<Real const>& cc_src_arr = cc_src.const_array(mfi);
149 const Array4<const Real> & stg_cons = S_stg_data[
IntVars::cons].const_array(mfi);
150 const Array4<const Real> & stg_xmom = S_stg_data[
IntVars::xmom].const_array(mfi);
151 const Array4<const Real> & stg_ymom = S_stg_data[
IntVars::ymom].const_array(mfi);
152 const Array4<const Real> & stg_zmom = S_stg_data[
IntVars::zmom].const_array(mfi);
153 const Array4<const Real> & prim = S_stg_prim.const_array(mfi);
154 const Array4<const Real> & qt_arr =
qt.const_array(mfi);
156 const Array4<const Real>& slow_rhs_cons = S_slow_rhs[
IntVars::cons].const_array(mfi);
157 const Array4<const Real>& slow_rhs_rho_u = S_slow_rhs[
IntVars::xmom].const_array(mfi);
158 const Array4<const Real>& slow_rhs_rho_v = S_slow_rhs[
IntVars::ymom].const_array(mfi);
159 const Array4<const Real>& slow_rhs_rho_w = S_slow_rhs[
IntVars::zmom].const_array(mfi);
161 const Array4<Real>& cur_cons = S_data[
IntVars::cons].array(mfi);
162 const Array4<Real>& cur_xmom = S_data[
IntVars::xmom].array(mfi);
163 const Array4<Real>& cur_ymom = S_data[
IntVars::ymom].array(mfi);
164 const Array4<Real>& cur_zmom = S_data[
IntVars::zmom].array(mfi);
166 const Array4<Real>& lagged = lagged_delta_rt.array(mfi);
168 const Array4<const Real>& prev_cons = S_prev[
IntVars::cons].const_array(mfi);
169 const Array4<const Real>& prev_xmom = S_prev[
IntVars::xmom].const_array(mfi);
170 const Array4<const Real>& prev_ymom = S_prev[
IntVars::ymom].const_array(mfi);
171 const Array4<const Real>& prev_zmom = S_prev[
IntVars::zmom].const_array(mfi);
174 const Array4<Real>& avg_xmom_arr = avg_xmom.array(mfi);
175 const Array4<Real>& avg_ymom_arr = avg_ymom.array(mfi);
176 const Array4<Real>& avg_zmom_arr = avg_zmom.array(mfi);
178 const Array4<const Real>& z_nd_old = z_phys_nd_old->const_array(mfi);
179 const Array4<const Real>& z_nd_new = z_phys_nd_new->const_array(mfi);
180 const Array4<const Real>& z_nd_stg = z_phys_nd_stg->const_array(mfi);
181 const Array4<const Real>& detJ_old = detJ_cc_old->const_array(mfi);
182 const Array4<const Real>& detJ_new = detJ_cc_new->const_array(mfi);
183 const Array4<const Real>& detJ_stg = detJ_cc_stg->const_array(mfi);
185 const Array4<const Real>& z_t_arr = z_t_rk->const_array(mfi);
186 const Array4<const Real>& zp_t_arr = z_t_pert->const_array(mfi);
188 const Array4< Real>& omega_arr = Omega.array(mfi);
190 const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
192 const Array4<Real>& theta_extrap = extrap.array(mfi);
195 const Array4<const Real>& mf_mx = mapfac[
MapFacType::m_x]->const_array(mfi);
196 const Array4<const Real>& mf_my = mapfac[
MapFacType::m_y]->const_array(mfi);
197 const Array4<const Real>& mf_ux = mapfac[
MapFacType::u_x]->const_array(mfi);
198 const Array4<const Real>& mf_vy = mapfac[
MapFacType::v_y]->const_array(mfi);
202 Box gbx = mfi.tilebox(); gbx.grow(1);
203 Box gtbx = mfi.nodaltilebox(0); gtbx.grow(1); gtbx.setSmall(2,0);
204 Box gtby = mfi.nodaltilebox(1); gtby.grow(1); gtby.setSmall(2,0);
208 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
213 theta_extrap(i,j,k) = delta_rt;
217 theta_extrap(i,j,k) *= (1.0 + RvOverRd*
qv);
220 lagged(i,j,k) = delta_rt;
222 }
else if (use_lagged_delta_rt) {
225 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
227 theta_extrap(i,j,k) = delta_rt + beta_d * (delta_rt - lagged(i,j,k));
231 theta_extrap(i,j,k) *= (1.0 + RvOverRd*
qv);
234 lagged(i,j,k) = delta_rt;
239 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
244 theta_extrap(i,j,k) *= (1.0 + RvOverRd*
qv);
248 RHS_fab.resize (tbz,1, The_Async_Arena());
249 soln_fab.resize (tbz,1, The_Async_Arena());
250 temp_rhs_fab.resize(tbz,2, The_Async_Arena());
252 auto const& RHS_a = RHS_fab.array();
253 auto const& soln_a = soln_fab.array();
254 auto const& temp_rhs_arr = temp_rhs_fab.array();
256 auto const& coeffA_a = coeff_A_mf.array(mfi);
257 auto const& inv_coeffB_a = inv_coeff_B_mf.array(mfi);
258 auto const& coeffC_a = coeff_C_mf.array(mfi);
259 auto const& coeffP_a = coeff_P_mf.array(mfi);
260 auto const& coeffQ_a = coeff_Q_mf.array(mfi);
266 BL_PROFILE(
"substep_xymom_T");
267 ParallelFor(tbx, tby,
268 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
273 Real gp_xi = (theta_extrap(i,j,k) - theta_extrap(i-1,j,k)) * dxi;
274 Real gp_zeta_on_iface = (k == 0) ?
275 0.5 * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
276 -theta_extrap(i-1,j,k ) - theta_extrap(i,j,k ) ) :
277 0.25 * dzi * ( theta_extrap(i-1,j,k+1) + theta_extrap(i,j,k+1)
278 -theta_extrap(i-1,j,k-1) - theta_extrap(i,j,k-1) );
279 Real gpx = h_zeta_old * gp_xi - h_xi_old * gp_zeta_on_iface;
282 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i-1,j,k) + qt_arr(i,j,k)) : 0.0;
284 Real pi_c = 0.5 * (pi_stage_ca(i-1,j,k,0) + pi_stage_ca(i ,j,k,0));
288 cur_xmom(i,j,k) = h_zeta_old * prev_xmom(i,j,k) + dtau * fast_rhs_rho_u
289 + dtau * slow_rhs_rho_u(i,j,k)
290 + dtau * xmom_src_arr(i,j,k);
292 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
297 Real gp_eta = (theta_extrap(i,j,k) -theta_extrap(i,j-1,k)) * dyi;
298 Real gp_zeta_on_jface = (k == 0) ?
299 0.5 * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
300 -theta_extrap(i,j,k ) - theta_extrap(i,j-1,k ) ) :
301 0.25 * dzi * ( theta_extrap(i,j,k+1) + theta_extrap(i,j-1,k+1)
302 -theta_extrap(i,j,k-1) - theta_extrap(i,j-1,k-1) );
303 Real gpy = h_zeta_old * gp_eta - h_eta_old * gp_zeta_on_jface;
306 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j-1,k) + qt_arr(i,j,k)) : 0.0;
308 Real pi_c = 0.5 * (pi_stage_ca(i,j-1,k,0) + pi_stage_ca(i,j ,k,0));
312 cur_ymom(i, j, k) = h_zeta_old * prev_ymom(i,j,k) + dtau * fast_rhs_rho_v
313 + dtau * slow_rhs_rho_v(i,j,k)
314 + dtau * ymom_src_arr(i,j,k);
321 for (
int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
322 flux[dir].resize(surroundingNodes(bx,dir),2,The_Async_Arena());
323 flux[dir].setVal<RunOn::Device>(0.);
325 const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
326 flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
330 BL_PROFILE(
"fast_T_making_rho_rhs");
331 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
335 Real xflux_lo = cur_xmom(i ,j,k) - stg_xmom(i ,j,k)*h_zeta_stg_xlo;
336 Real xflux_hi = cur_xmom(i+1,j,k) - stg_xmom(i+1,j,k)*h_zeta_stg_xhi;
340 Real yflux_lo = cur_ymom(i,j ,k) - stg_ymom(i,j ,k)*h_zeta_stg_ylo;
341 Real yflux_hi = cur_ymom(i,j+1,k) - stg_ymom(i,j+1,k)*h_zeta_stg_yhi;
344 temp_rhs_arr(i,j,k,0) = ( ( xflux_hi - xflux_lo ) * dxi + ( yflux_hi - yflux_lo ) * dyi );
345 temp_rhs_arr(i,j,k,1) = ( (( xflux_hi * (prim(i,j,k,0) + prim(i+1,j,k,0)) -
346 xflux_lo * (prim(i,j,k,0) + prim(i-1,j,k,0)) ) * dxi +
347 ( yflux_hi * (prim(i,j,k,0) + prim(i,j+1,k,0)) -
348 yflux_lo * (prim(i,j,k,0) + prim(i,j-1,k,0)) ) * dyi) * 0.5 );
351 (flx_arr[0])(i,j,k,0) = xflux_lo;
352 (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));
354 (flx_arr[1])(i,j,k,0) = yflux_lo;
355 (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));
358 (flx_arr[0])(i+1,j,k,0) = xflux_hi;
359 (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));
362 (flx_arr[1])(i,j+1,k,0) = yflux_hi;
363 (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));
373 Box gbxo = mfi.nodaltilebox(2);
375 BL_PROFILE(
"fast_MT_making_omega");
376 Box gbxo_lo = gbxo; gbxo_lo.setBig(2,0);
377 ParallelFor(gbxo_lo, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
378 omega_arr(i,j,k) = 0.;
380 Box gbxo_hi = gbxo; gbxo_hi.setSmall(2,gbxo.bigEnd(2));
381 ParallelFor(gbxo_hi, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
382 omega_arr(i,j,k) = prev_zmom(i,j,k) - stg_zmom(i,j,k) - zp_t_arr(i,j,k);
384 Box gbxo_mid = gbxo; gbxo_mid.setSmall(2,1); gbxo_mid.setBig(2,gbxo.bigEnd(2)-1);
385 ParallelFor(gbxo_mid, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
387 (
OmegaFromW(i,j,k,prev_zmom(i,j,k),prev_xmom,prev_ymom,mf_ux,mf_vy,z_nd_old,dxInv)
388 -
OmegaFromW(i,j,k, stg_zmom(i,j,k), stg_xmom, stg_ymom,mf_ux,mf_vy,z_nd_old,dxInv) )
394 ParallelFor(tbx, tby,
395 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
398 cur_xmom(i, j, k) /= h_zeta_new;
399 avg_xmom_arr(i,j,k) += facinv*(cur_xmom(i,j,k) - stg_xmom(i,j,k));
401 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
404 cur_ymom(i, j, k) /= h_zeta_new;
405 avg_ymom_arr(i,j,k) += facinv*(cur_ymom(i,j,k) - stg_ymom(i,j,k));
408 Box bx_shrunk_in_k = bx;
409 int klo = tbz.smallEnd(2);
410 int khi = tbz.bigEnd(2);
411 bx_shrunk_in_k.setSmall(2,klo+1);
412 bx_shrunk_in_k.setBig(2,khi-1);
417 Real halfg = std::abs(0.5 * grav_gpu[2]);
420 BL_PROFILE(
"fast_loop_on_shrunk_t");
422 ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
424 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k-1) + qt_arr(i,j,k)) : 0.0;
426 Real coeff_P = coeffP_a(i,j,k) / (1.0 + q);
427 Real coeff_Q = coeffQ_a(i,j,k) / (1.0 + q);
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);
451 Real detJdiff = (detJ_old(i,j,k) - detJ_old(i,j,k-1)) / (detJ_old(i,j,k)*detJ_old(i,j,k-1));
454 R1_tmp += halfg * ( beta_1 * dzi * (Omega_kp1/detJ_old(i,j,k) + detJdiff*Omega_k - Omega_km1/detJ_old(i,j,k-1))
455 + 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) );
458 R1_tmp += -( coeff_P/detJ_stg(i,j,k ) * ( beta_1 * dzi * (Omega_kp1*theta_t_hi - Omega_k*theta_t_mid) +
460 + coeff_Q/detJ_stg(i,j,k-1) * ( beta_1 * dzi * (Omega_k*theta_t_mid - Omega_km1*theta_t_lo) +
464 Real detJ_half = 0.5 * (detJ_stg(i,j,k) + detJ_stg(i,j,k-1));
465 RHS_a(i,j,k) = prev_zmom(i,j,k) - stg_zmom(i,j,k)
466 + dtau * (slow_rhs_rho_w(i,j,k) + zmom_src_arr(i,j,k)) / detJ_half
467 + dtau * (R0_tmp + dtau*beta_2*R1_tmp);
470 Real UppVpp =
OmegaFromW(i,j,k,0.,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
471 -
OmegaFromW(i,j,k,0.,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
472 RHS_a(i,j,k) += UppVpp;
479 auto const lo = lbound(bx);
480 auto const hi = ubound(bx);
483 BL_PROFILE(
"substep_b2d_loop_t");
486 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
489 Real rho_on_bdy = 0.5 * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
490 RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,0);
492 soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
495 RHS_a(i,j,hi.z+1) = 0.0;
497 for (
int k = lo.z+1; k <= hi.z+1; k++) {
498 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);
501 for (
int k = hi.z; k >= lo.z; k--) {
502 soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
506 cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
509 for (
int j = lo.y; j <= hi.y; ++j) {
511 for (
int i = lo.x; i <= hi.x; ++i) {
513 Real rho_on_bdy = 0.5 * ( prev_cons(i,j,lo.z) + prev_cons(i,j,lo.z-1) );
514 RHS_a(i,j,lo.z) = rho_on_bdy * zp_t_arr(i,j,lo.z);
516 soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
520 for (
int j = lo.y; j <= hi.y; ++j) {
522 for (
int i = lo.x; i <= hi.x; ++i) {
523 RHS_a (i,j,hi.z+1) = 0.0;
526 for (
int k = lo.z+1; k <= hi.z+1; ++k) {
527 for (
int j = lo.y; j <= hi.y; ++j) {
529 for (
int i = lo.x; i <= hi.x; ++i) {
530 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);
534 for (
int k = hi.z; k >= lo.z; --k) {
535 for (
int j = lo.y; j <= hi.y; ++j) {
537 for (
int i = lo.x; i <= hi.x; ++i) {
538 soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) * soln_a(i,j,k+1);
544 for (
int j = lo.y; j <= hi.y; ++j) {
546 for (
int i = lo.x; i <= hi.x; ++i) {
547 cur_zmom(i,j,hi.z+1) = stg_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
554 BL_PROFILE(
"substep_new_drhow");
556 ParallelFor(tbz, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
561 cur_zmom(i,j,k) =
WFromOmega(i,j,k,rho_on_face*(z_t_arr(i,j,k)+zp_t_arr(i,j,k)),
562 cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv);
569 Real UppVpp =
WFromOmega(i,j,k,0.0,cur_xmom,cur_ymom,mf_ux,mf_vy,z_nd_new,dxInv)
570 -
WFromOmega(i,j,k,0.0,stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
571 Real wpp = soln_a(i,j,k) + UppVpp;
572 Real dJ_old_kface = 0.5 * (detJ_old(i,j,k) + detJ_old(i,j,k-1));
573 Real dJ_new_kface = 0.5 * (detJ_new(i,j,k) + detJ_new(i,j,k-1));
575 cur_zmom(i,j,k) = dJ_old_kface * (stg_zmom(i,j,k) + wpp);
576 cur_zmom(i,j,k) /= dJ_new_kface;
578 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)
579 -
OmegaFromW(i,j,k,stg_zmom(i,j,k),stg_xmom,stg_ymom,mf_ux,mf_vy,z_nd_stg,dxInv);
580 soln_a(i,j,k) -= rho_on_face * zp_t_arr(i,j,k);
589 BL_PROFILE(
"fast_rho_final_update");
590 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
592 Real zflux_lo = beta_2 * soln_a(i,j,k ) + beta_1 * omega_arr(i,j,k);
593 Real zflux_hi = beta_2 * soln_a(i,j,k+1) + beta_1 * omega_arr(i,j,k+1);
597 avg_zmom_arr(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
599 (flx_arr[2])(i,j,k,0) = zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
605 Real fast_rhs_rho = -(temp_rhs_arr(i,j,k,0) + ( zflux_hi - zflux_lo ) * dzi);
606 Real temp_rho = detJ_old(i,j,k) * cur_cons(i,j,k,0) +
607 dtau * ( slow_rhs_cons(i,j,k,0)*detJ_stg(i,j,k) + fast_rhs_rho );
608 cur_cons(i,j,k,0) = temp_rho / detJ_new(i,j,k);
613 Real fast_rhs_rhotheta = -( temp_rhs_arr(i,j,k,1) + 0.5 *
614 ( zflux_hi * (prim(i,j,k) + prim(i,j,k+1))
615 - zflux_lo * (prim(i,j,k) + prim(i,j,k-1)) ) * dzi );
616 Real temp_rth = detJ_old(i,j,k) * cur_cons(i,j,k,1) +
617 dtau * ( slow_rhs_cons(i,j,k,1)*detJ_stg(i,j,k) + fast_rhs_rhotheta );
618 cur_cons(i,j,k,1) = temp_rth / detJ_new(i,j,k);
620 (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));
624 avg_zmom_arr(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
626 (flx_arr[2])(i,j,k+1,0) = zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
627 (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));
639 int strt_comp_reflux = 0;
641 int num_comp_reflux = 1;
642 if (level < finest_level) {
643 fr_as_crse->CrseAdd(mfi,
644 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
645 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
648 fr_as_fine->FineAdd(mfi,
649 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
650 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
656 Gpu::streamSynchronize();
constexpr amrex::Real R_v
Definition: ERF_Constants.H:11
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:24
@ m_y
Definition: ERF_DataStruct.H:24
@ u_x
Definition: ERF_DataStruct.H:23
@ m_x
Definition: ERF_DataStruct.H:23
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:53
#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::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:412
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:115
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:102
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:142
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:168
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:462
@ 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
@ qt
Definition: ERF_Kessler.H:27
@ qv
Definition: ERF_Kessler.H:28