74 BL_PROFILE_REGION(
"erf_fast_rhs_N()");
76 Real beta_1 = 0.5 * (1.0 - beta_s);
77 Real beta_2 = 0.5 * (1.0 + beta_s);
84 const Real* dx = geom.CellSize();
85 const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
92 const auto& dm = S_stage_data[
IntVars::cons].DistributionMap();
94 MultiFab Delta_rho_theta( ba , dm, 1, 1);
95 MultiFab Delta_rho_w (convert(ba,IntVect(0,0,1)), dm, 1, IntVect(1,1,0));
97 MultiFab coeff_A_mf(fast_coeffs, make_alias, 0, 1);
98 MultiFab inv_coeff_B_mf(fast_coeffs, make_alias, 1, 1);
99 MultiFab coeff_C_mf(fast_coeffs, make_alias, 2, 1);
100 MultiFab coeff_P_mf(fast_coeffs, make_alias, 3, 1);
101 MultiFab coeff_Q_mf(fast_coeffs, make_alias, 4, 1);
105 const Array<Real,AMREX_SPACEDIM> grav{0.0, 0.0, -gravity};
106 const GpuArray<Real,AMREX_SPACEDIM> grav_gpu{grav[0], grav[1], grav[2]};
119 AMREX_ALWAYS_ASSERT(nrk > 0 || step == 0);
126 #pragma omp parallel if (Gpu::notInLaunchRegion())
128 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
130 const Array4<const Real>& prev_cons = S_prev[
IntVars::cons].const_array(mfi);
131 const Array4<const Real>& prev_zmom = S_prev[
IntVars::zmom].const_array(mfi);
133 const Array4<const Real>& stage_cons = S_stage_data[
IntVars::cons].const_array(mfi);
134 const Array4<const Real>& stage_zmom = S_stage_data[
IntVars::zmom].const_array(mfi);
136 const Array4<Real>& prev_drho_w = Delta_rho_w.array(mfi);
137 const Array4<Real>& prev_drho_theta = Delta_rho_theta.array(mfi);
138 const Array4<Real>& lagged_arr = lagged_delta_rt.array(mfi);
139 const Array4<Real>& theta_extrap = extrap.array(mfi);
140 const Array4<const Real>& prim = S_stage_prim.const_array(mfi);
142 Box gbx = mfi.growntilebox(1);
143 ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
148 theta_extrap(i,j,k) = prev_drho_theta(i,j,k);
150 theta_extrap(i,j,k) = prev_drho_theta(i,j,k) + beta_d *
151 ( prev_drho_theta(i,j,k) - lagged_arr(i,j,k) );
156 theta_extrap(i,j,k) *= (1.0 + RvOverRd*
qv);
160 lagged_arr(i,j,k) = prev_drho_theta(i,j,k);
165 Box gtbz = mfi.nodaltilebox(2);
166 gtbz.grow(IntVect(1,1,0));
167 ParallelFor(gtbz, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
168 prev_drho_w(i,j,k) = prev_zmom(i,j,k) - stage_zmom(i,j,k);
177 #pragma omp parallel if (Gpu::notInLaunchRegion())
179 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
181 Box tbx = mfi.nodaltilebox(0);
182 Box tby = mfi.nodaltilebox(1);
184 const Array4<Real const>& xmom_src_arr = xmom_src.const_array(mfi);
185 const Array4<Real const>& ymom_src_arr = ymom_src.const_array(mfi);
187 const Array4<const Real> & stage_xmom = S_stage_data[
IntVars::xmom].const_array(mfi);
188 const Array4<const Real> & stage_ymom = S_stage_data[
IntVars::ymom].const_array(mfi);
189 const Array4<const Real> & qt_arr =
qt.const_array(mfi);
191 const Array4<const Real>& slow_rhs_rho_u = S_slow_rhs[
IntVars::xmom].const_array(mfi);
192 const Array4<const Real>& slow_rhs_rho_v = S_slow_rhs[
IntVars::ymom].const_array(mfi);
194 const Array4<Real>& temp_cur_xmom_arr = temp_cur_xmom.array(mfi);
195 const Array4<Real>& temp_cur_ymom_arr = temp_cur_ymom.array(mfi);
197 const Array4<const Real>& prev_xmom = S_prev[
IntVars::xmom].const_array(mfi);
198 const Array4<const Real>& prev_ymom = S_prev[
IntVars::ymom].const_array(mfi);
201 const Array4< Real>& avg_xmom_arr = avg_xmom.array(mfi);
202 const Array4< Real>& avg_ymom_arr = avg_ymom.array(mfi);
204 const Array4<const Real>& pi_stage_ca = pi_stage.const_array(mfi);
206 const Array4<Real>& theta_extrap = extrap.array(mfi);
209 const Array4<const Real>& mf_ux = mapfac[
MapFacType::u_x]->const_array(mfi);
210 const Array4<const Real>& mf_vy = mapfac[
MapFacType::v_y]->const_array(mfi);
215 if (nrk == 0 and step == 0) {
216 ParallelFor(tbx, tby,
217 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
219 Real new_drho_u = dtau * slow_rhs_rho_u(i,j,k) + dtau * xmom_src_arr(i,j,k);;
220 avg_xmom_arr(i,j,k) += facinv*new_drho_u;
221 temp_cur_xmom_arr(i,j,k) = stage_xmom(i,j,k) + new_drho_u;
223 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
225 Real new_drho_v = dtau * slow_rhs_rho_v(i,j,k) + dtau * ymom_src_arr(i,j,k);
226 avg_ymom_arr(i,j,k) += facinv*new_drho_v;
227 temp_cur_ymom_arr(i,j,k) = stage_ymom(i,j,k) + new_drho_v;
230 ParallelFor(tbx, tby,
231 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
234 Real gpx = (theta_extrap(i,j,k) - theta_extrap(i-1,j,k))*dxi;
237 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i-1,j,k)) : 0.0;
239 Real pi_c = 0.5 * (pi_stage_ca(i-1,j,k,0) + pi_stage_ca(i,j,k,0));
242 Real new_drho_u = prev_xmom(i,j,k) - stage_xmom(i,j,k)
243 + dtau * fast_rhs_rho_u + dtau * slow_rhs_rho_u(i,j,k)
244 + dtau * xmom_src_arr(i,j,k);
246 avg_xmom_arr(i,j,k) += facinv*new_drho_u;
248 temp_cur_xmom_arr(i,j,k) = stage_xmom(i,j,k) + new_drho_u;
250 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
253 Real gpy = (theta_extrap(i,j,k) - theta_extrap(i,j-1,k))*dyi;
256 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i,j-1,k)) : 0.0;
258 Real pi_c = 0.5 * (pi_stage_ca(i,j-1,k,0) + pi_stage_ca(i,j,k,0));
261 Real new_drho_v = prev_ymom(i,j,k) - stage_ymom(i,j,k)
262 + dtau * fast_rhs_rho_v + dtau * slow_rhs_rho_v(i,j,k)
263 + dtau * ymom_src_arr(i,j,k);
265 avg_ymom_arr(i,j,k) += facinv*new_drho_v;
267 temp_cur_ymom_arr(i,j,k) = stage_ymom(i,j,k) + new_drho_v;
273 #pragma omp parallel if (Gpu::notInLaunchRegion())
276 std::array<FArrayBox,AMREX_SPACEDIM> flux;
279 Box bx = mfi.tilebox();
280 Box tbz = surroundingNodes(bx,2);
282 Box vbx = mfi.validbox();
283 const auto& vbx_hi = ubound(vbx);
285 const Array4<Real const>& zmom_src_arr = zmom_src.const_array(mfi);
287 const Array4<const Real>& stage_xmom = S_stage_data[
IntVars::xmom].const_array(mfi);
288 const Array4<const Real>& stage_ymom = S_stage_data[
IntVars::ymom].const_array(mfi);
289 const Array4<const Real>& stage_zmom = S_stage_data[
IntVars::zmom].const_array(mfi);
290 const Array4<const Real> & prim = S_stage_prim.const_array(mfi);
291 const Array4<const Real> & qt_arr =
qt.const_array(mfi);
293 const Array4<const Real>& prev_drho_theta = Delta_rho_theta.array(mfi);
295 const Array4<const Real>& prev_cons = S_prev[
IntVars::cons].const_array(mfi);
296 const Array4<const Real>& stage_cons = S_stage_data[
IntVars::cons].const_array(mfi);
298 const Array4<const Real>& slow_rhs_cons = S_slow_rhs[
IntVars::cons].const_array(mfi);
299 const Array4<const Real>& slow_rhs_rho_w = S_slow_rhs[
IntVars::zmom].const_array(mfi);
301 const Array4<const Real>& prev_zmom = S_prev[
IntVars::zmom].const_array(mfi);
302 const Array4< Real>& cur_zmom = S_data[
IntVars::zmom].array(mfi);
304 const Array4<Real>& temp_cur_xmom_arr = temp_cur_xmom.array(mfi);
305 const Array4<Real>& temp_cur_ymom_arr = temp_cur_ymom.array(mfi);
308 const Array4< Real>& avg_zmom_arr = avg_zmom.array(mfi);
311 const Array4<const Real>& mf_mx = mapfac[
MapFacType::m_x]->const_array(mfi);
312 const Array4<const Real>& mf_my = mapfac[
MapFacType::m_y]->const_array(mfi);
313 const Array4<const Real>& mf_uy = mapfac[
MapFacType::u_y]->const_array(mfi);
314 const Array4<const Real>& mf_vx = mapfac[
MapFacType::v_x]->const_array(mfi);
317 RHS_fab.resize(tbz,1, The_Async_Arena());
320 soln_fab.resize(tbz,1, The_Async_Arena());
322 auto const& RHS_a = RHS_fab.array();
323 auto const& soln_a = soln_fab.array();
325 auto const& temp_rhs_arr = temp_rhs.array(mfi);
327 auto const& coeffA_a = coeff_A_mf.array(mfi);
328 auto const& inv_coeffB_a = inv_coeff_B_mf.array(mfi);
329 auto const& coeffC_a = coeff_C_mf.array(mfi);
330 auto const& coeffP_a = coeff_P_mf.array(mfi);
331 auto const& coeffQ_a = coeff_Q_mf.array(mfi);
336 for (
int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
337 flux[dir].resize(surroundingNodes(bx,dir),2);
338 flux[dir].setVal<RunOn::Device>(0.);
340 const GpuArray<const Array4<Real>, AMREX_SPACEDIM>
341 flx_arr{{AMREX_D_DECL(flux[0].array(), flux[1].array(), flux[2].array())}};
344 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
345 Real xflux_lo = (temp_cur_xmom_arr(i ,j,k) - stage_xmom(i ,j,k)) / mf_uy(i ,j,0);
346 Real xflux_hi = (temp_cur_xmom_arr(i+1,j,k) - stage_xmom(i+1,j,k)) / mf_uy(i+1,j,0);
347 Real yflux_lo = (temp_cur_ymom_arr(i,j ,k) - stage_ymom(i,j ,k)) / mf_vx(i,j ,0);
348 Real yflux_hi = (temp_cur_ymom_arr(i,j+1,k) - stage_ymom(i,j+1,k)) / mf_vx(i,j+1,0);
350 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
352 temp_rhs_arr(i,j,k,
Rho_comp ) = ( xflux_hi - xflux_lo ) * dxi * mfsq
353 + ( yflux_hi - yflux_lo ) * dyi * mfsq;
354 temp_rhs_arr(i,j,k,
RhoTheta_comp) = (( xflux_hi * (prim(i,j,k,0) + prim(i+1,j,k,0)) -
355 xflux_lo * (prim(i,j,k,0) + prim(i-1,j,k,0)) ) * dxi * mfsq +
356 ( yflux_hi * (prim(i,j,k,0) + prim(i,j+1,k,0)) -
357 yflux_lo * (prim(i,j,k,0) + prim(i,j-1,k,0)) ) * dyi * mfsq) * 0.5;
360 (flx_arr[0])(i,j,k,0) = xflux_lo;
361 (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));
363 (flx_arr[1])(i,j,k,0) = yflux_lo;
364 (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));
367 (flx_arr[0])(i+1,j,k,0) = xflux_hi;
368 (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));
371 (flx_arr[1])(i,j+1,k,0) = yflux_hi;
372 (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));
377 Box bx_shrunk_in_k = bx;
378 int klo = tbz.smallEnd(2);
379 int khi = tbz.bigEnd(2);
380 bx_shrunk_in_k.setSmall(2,klo+1);
381 bx_shrunk_in_k.setBig(2,khi-1);
386 Real halfg = std::abs(0.5 * grav_gpu[2]);
392 ParallelFor(bx_shrunk_in_k, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
394 Real q = (l_use_moisture) ? 0.5 * (qt_arr(i,j,k) + qt_arr(i,j,k-1)) : 0.0;
396 Real coeff_P = coeffP_a(i,j,k) / (1.0 + q);
397 Real coeff_Q = coeffQ_a(i,j,k) / (1.0 + q);
403 Real Omega_kp1 = prev_zmom(i,j,k+1) - stage_zmom(i,j,k+1);
404 Real Omega_k = prev_zmom(i,j,k ) - stage_zmom(i,j,k );
405 Real Omega_km1 = prev_zmom(i,j,k-1) - stage_zmom(i,j,k-1);
410 Real R0_tmp = coeff_P * prev_drho_theta(i,j,k) + coeff_Q * prev_drho_theta(i,j,k-1)
411 - halfg * ( old_drho_k + old_drho_km1 );
416 +temp_rhs_arr(i,j,k,0) + temp_rhs_arr(i,j,k-1) )
421 R1_tmp += beta_1 * dzi * ( (Omega_kp1 - Omega_km1) * halfg
422 -(Omega_kp1*theta_t_hi - Omega_k *theta_t_mid) * coeff_P
423 -(Omega_k *theta_t_mid - Omega_km1*theta_t_lo ) * coeff_Q );
426 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));
433 auto const lo = lbound(bx);
434 auto const hi = ubound(bx);
436 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
439 RHS_a (i,j,lo.z) = prev_zmom(i,j,lo.z) - stage_zmom(i,j,lo.z)
440 + dtau * slow_rhs_rho_w(i,j,lo.z)
441 + dtau * zmom_src_arr(i,j,lo.z);
444 RHS_a (i,j,hi.z+1) = prev_zmom(i,j,hi.z+1) - stage_zmom(i,j,hi.z+1)
445 + dtau * slow_rhs_rho_w(i,j,hi.z+1)
446 + dtau * zmom_src_arr(i,j,hi.z+1);
450 if (l_implicit_substepping) {
452 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
455 soln_a(i,j,lo.z) = RHS_a(i,j,lo.z) * inv_coeffB_a(i,j,lo.z);
456 cur_zmom(i,j,lo.z) = stage_zmom(i,j,lo.z) + soln_a(i,j,lo.z);
458 for (
int k = lo.z+1; k <= hi.z+1; k++) {
459 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);
462 cur_zmom(i,j,hi.z+1) = stage_zmom(i,j,hi.z+1) + soln_a(i,j,hi.z+1);
464 for (
int k = hi.z; k >= lo.z; k--) {
465 soln_a(i,j,k) -= ( coeffC_a(i,j,k) * inv_coeffB_a(i,j,k) ) *soln_a(i,j,k+1);
466 cur_zmom(i,j,k) = stage_zmom(i,j,k) + soln_a(i,j,k);
472 ParallelFor(b2d, [=] AMREX_GPU_DEVICE (
int i,
int j,
int)
474 for (
int k = lo.z; k <= hi.z+1; k++) {
475 soln_a(i,j,k) = RHS_a(i,j,k);
476 cur_zmom(i,j,k) = stage_zmom(i,j,k) + soln_a(i,j,k);
481 if (l_implicit_substepping) {
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);
514 for (
int k = lo.z; k <= hi.z+1; ++k) {
515 for (
int j = lo.y; j <= hi.y; ++j) {
517 for (
int i = lo.x; i <= hi.x; ++i) {
519 soln_a(i,j,k) = RHS_a(i,j,k);
521 cur_zmom(i,j,k) = stage_zmom(i,j,k) + soln_a(i,j,k);
532 const Array4<Real>& prev_drho_w = Delta_rho_w.array(mfi);
533 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
535 Real zflux_lo = beta_2 * soln_a(i,j,k ) + beta_1 * prev_drho_w(i,j,k );
536 Real zflux_hi = beta_2 * soln_a(i,j,k+1) + beta_1 * prev_drho_w(i,j,k+1);
538 avg_zmom_arr(i,j,k) += facinv*zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
540 (flx_arr[2])(i,j,k,0) = zflux_lo / (mf_mx(i,j,0) * mf_my(i,j,0));
541 (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));
545 avg_zmom_arr(i,j,k+1) += facinv * zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
547 (flx_arr[2])(i,j,k+1,0) = zflux_hi / (mf_mx(i,j,0) * mf_my(i,j,0));
548 (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));
552 temp_rhs_arr(i,j,k,
Rho_comp ) += dzi * ( zflux_hi - zflux_lo );
553 temp_rhs_arr(i,j,k,
RhoTheta_comp) += 0.5 * dzi * ( zflux_hi * (prim(i,j,k) + prim(i,j,k+1))
554 - zflux_lo * (prim(i,j,k) + prim(i,j,k-1)) );
559 int strt_comp_reflux = 0;
561 int num_comp_reflux = 1;
562 if (level < finest_level) {
563 fr_as_crse->CrseAdd(mfi,
564 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
565 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
568 fr_as_fine->FineAdd(mfi,
569 {{AMREX_D_DECL(&(flux[0]), &(flux[1]), &(flux[2]))}},
570 dx, dtau, strt_comp_reflux, strt_comp_reflux, num_comp_reflux, RunOn::Device);
576 Gpu::streamSynchronize();
583 #pragma omp parallel if (Gpu::notInLaunchRegion())
585 for ( MFIter mfi(S_stage_data[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
587 const Box& bx = mfi.tilebox();
589 const Array4< Real>& cur_cons = S_data[
IntVars::cons].array(mfi);
590 const Array4<const Real>& prev_cons = S_prev[
IntVars::cons].const_array(mfi);
591 auto const& temp_rhs_arr = temp_rhs.const_array(mfi);
592 auto const& slow_rhs_cons = S_slow_rhs[
IntVars::cons].const_array(mfi);
593 const Array4<Real const>& cc_src_arr = cc_src.const_array(mfi);
596 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
608 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
622 const Array4<Real>& cur_xmom = S_data[
IntVars::xmom].array(mfi);
623 const Array4<Real>& cur_ymom = S_data[
IntVars::ymom].array(mfi);
625 const Array4<Real const>& temp_cur_xmom_arr = temp_cur_xmom.const_array(mfi);
626 const Array4<Real const>& temp_cur_ymom_arr = temp_cur_ymom.const_array(mfi);
628 Box tbx = surroundingNodes(bx,0);
629 Box tby = surroundingNodes(bx,1);
631 ParallelFor(tbx, tby,
632 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
634 cur_xmom(i,j,k) = temp_cur_xmom_arr(i,j,k);
636 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
638 cur_ymom(i,j,k) = temp_cur_ymom_arr(i,j,k);
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_x
Definition: ERF_DataStruct.H:22
@ u_y
Definition: ERF_DataStruct.H:23
@ v_y
Definition: ERF_DataStruct.H:23
@ m_y
Definition: ERF_DataStruct.H:23
@ u_x
Definition: ERF_DataStruct.H:22
@ m_x
Definition: ERF_DataStruct.H:22
#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_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
@ 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