7 template<
typename InterpType_H,
typename InterpType_V,
typename WallInterpType>
10 const amrex::Vector<amrex::Box>& bxy_grown,
11 const amrex::Vector<amrex::Box>& bxz_grown,
12 const amrex::Array4<const amrex::Real>& rho_u,
13 const amrex::Array4<const amrex::Real>& rho_v,
14 const amrex::Array4<const amrex::Real>& rho_w,
15 const amrex::Array4<const amrex::Real>& u,
16 const amrex::Array4<const amrex::Real>& v,
17 const amrex::Array4<const amrex::Real>& w,
18 const amrex::Array4<const amrex::EBCellFlag>& u_cflag,
19 const amrex::Array4<const amrex::Real>& u_afrac_x,
20 const amrex::Array4<const amrex::Real>& u_afrac_y,
21 const amrex::Array4<const amrex::Real>& u_afrac_z,
22 const amrex::Array4<const amrex::EBCellFlag>& v_cflag,
23 const amrex::Array4<const amrex::Real>& v_afrac_x,
24 const amrex::Array4<const amrex::Real>& v_afrac_y,
25 const amrex::Array4<const amrex::Real>& v_afrac_z,
26 const amrex::Array4<const amrex::EBCellFlag>& w_cflag,
27 const amrex::Array4<const amrex::Real>& w_afrac_x,
28 const amrex::Array4<const amrex::Real>& w_afrac_y,
29 const amrex::Array4<const amrex::Real>& w_afrac_z,
30 const amrex::Array4<const amrex::Real>& mf_ux_inv,
31 const amrex::Array4<const amrex::Real>& mf_vx_inv,
32 const amrex::Array4<const amrex::Real>& mf_uy_inv,
33 const amrex::Array4<const amrex::Real>& mf_vy_inv,
34 const amrex::Real upw_frac_h,
35 const amrex::Real upw_frac_v,
37 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_u_arr,
38 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_v_arr,
39 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_w_arr,
40 const int lo_z_face,
const int hi_z_face)
43 InterpType_H interp_u_horz(u, upw_frac_h); InterpType_V interp_u_vert(u, upw_frac_v);
44 InterpType_H interp_v_horz(v, upw_frac_h); InterpType_V interp_v_vert(v, upw_frac_v);
45 InterpType_H interp_w_horz(w, upw_frac_h); InterpType_V interp_w_vert(w, upw_frac_v);
48 const int ncells_u_horz = interp_u_horz.GetUpwindCellNumber();
49 const int ncells_u_vert = interp_u_vert.GetUpwindCellNumber();
50 const int ncells_v_horz = interp_v_horz.GetUpwindCellNumber();
51 const int ncells_v_vert = interp_v_vert.GetUpwindCellNumber();
52 const int ncells_w_horz = interp_w_horz.GetUpwindCellNumber();
53 const int ncells_w_vert = interp_w_vert.GetUpwindCellNumber();
60 UPWIND3 interp_u_UPW3(u, upw_frac_h);
61 UPWIND3 interp_v_UPW3(v, upw_frac_h);
62 UPWIND3 interp_w_UPW3(w, upw_frac_h);
68 amrex::ParallelFor(bxx_grown[0], bxx_grown[1], bxx_grown[2],
70 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
72 if ( u_afrac_x(i,j,k) > 0. )
74 amrex::Real rho_u_avg = 0.5 * (rho_u(i-1,j,k) * mf_ux_inv(i-1,j,0) + rho_u(i,j,k) * mf_ux_inv(i,j,0));
80 for (
int ii=0; ii<ncells_u_horz; ++ii) {
81 if (u_cflag(i-ii-1,j,k).isCovered()) {
87 else if (rho_u_avg<0.)
89 for (
int ii=0; ii<ncells_u_horz; ++ii) {
90 if (u_cflag(i+ii,j,k).isCovered()) {
99 amrex::Real interpx = 0.;
100 if (icell==ncells_u_horz) {
101 interp_u_horz.InterpolateInX(i,j,k,0,interpx,rho_u_avg);
105 }
else if (icell==2) {
109 flx_u_arr[0](i,j,k) = rho_u_avg * interpx;
113 flx_u_arr[0](i,j,k) = 0.;
117 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
119 if ( u_afrac_y(i,j,k) > 0. )
121 amrex::Real rho_v_avg = 0.5 * (rho_v(i,j,k) * mf_vx_inv(i,j,0) + rho_v(i-1,j,k) * mf_vx_inv(i-1,j,0));
127 for (
int jj=0; jj<ncells_u_horz; ++jj) {
128 if (u_cflag(i,j-jj-1,k).isCovered()) {
134 else if (rho_v_avg<0.)
136 for (
int jj=0; jj<ncells_u_horz; ++jj) {
137 if (u_cflag(i,j+jj,k).isCovered()) {
146 amrex::Real interpy = 0.;
147 if (jcell==ncells_u_horz) {
148 interp_u_horz.InterpolateInY(i,j,k,0,interpy,rho_v_avg);
152 }
else if (jcell==2) {
156 flx_u_arr[1](i,j,k) = rho_v_avg * interpy;
160 flx_u_arr[1](i,j,k) = 0.;
164 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
166 if ( u_afrac_z(i,j,k) > 0. )
168 amrex::Real rho_w_avg = 0.5 * (rho_w(i,j,k) + rho_w(i-1,j,k));
174 for (
int kk=0; kk<ncells_u_vert; ++kk) {
175 if (u_cflag(i,j,k-kk-1).isCovered()) {
181 else if (rho_w_avg<0.)
183 for (
int kk=0; kk<ncells_u_vert; ++kk) {
184 if (u_cflag(i,j,k+kk).isCovered()) {
193 amrex::Real interpz = 0.;
194 if (kcell==ncells_u_vert) {
195 interp_u_vert.InterpolateInZ(i,j,k,0,interpz,rho_w_avg);
199 }
else if (kcell==2) {
203 flx_u_arr[2](i,j,k) = rho_w_avg * interpz;
207 flx_u_arr[2](i,j,k) = 0.;
215 amrex::ParallelFor(bxy_grown[0], bxy_grown[1], bxy_grown[2],
217 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
219 if ( v_afrac_x(i,j,k) > 0. )
221 amrex::Real rho_u_avg = 0.5 * (rho_u(i,j,k) * mf_uy_inv(i,j,0) + rho_u(i,j-1,k) * mf_uy_inv(i,j-1,0));
227 for (
int ii=0; ii<ncells_v_horz; ++ii) {
228 if (v_cflag(i-ii-1,j,k).isCovered()) {
234 else if (rho_u_avg<0.)
236 for (
int ii=0; ii<ncells_v_horz; ++ii) {
237 if (v_cflag(i+ii,j,k).isCovered()) {
246 amrex::Real interpx = 0.;
247 if (icell==ncells_v_horz) {
248 interp_v_horz.InterpolateInX(i,j,k,0,interpx,rho_u_avg);
252 }
else if (icell==2) {
256 flx_v_arr[0](i,j,k) = rho_u_avg * interpx;
260 flx_v_arr[0](i,j,k) = 0.;
264 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
266 if ( v_afrac_y(i,j,k) > 0. )
268 amrex::Real rho_v_avg = 0.5 * (rho_v(i,j,k) * mf_vy_inv(i,j,0) + rho_v(i,j-1,k) * mf_vy_inv(i,j-1,0));
274 for (
int jj=0; jj<ncells_v_horz; ++jj) {
275 if (v_cflag(i,j-jj-1,k).isCovered()) {
280 else if (rho_v_avg<0.)
282 for (
int jj=0; jj<ncells_v_horz; ++jj) {
283 if (v_cflag(i,j+jj,k).isCovered()) {
291 amrex::Real interpy = 0.;
292 if (jcell==ncells_v_horz) {
293 interp_v_horz.InterpolateInY(i,j,k,0,interpy,rho_v_avg);
297 }
else if (jcell==2) {
301 flx_v_arr[1](i,j,k) = rho_v_avg * interpy;
305 flx_v_arr[1](i,j,k) = 0.;
309 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
311 if ( v_afrac_z(i,j,k) > 0. )
313 amrex::Real rho_w_avg = 0.5 * (rho_w(i,j,k) + rho_w(i,j-1,k));
319 for (
int kk=0; kk<ncells_v_vert; ++kk) {
320 if (v_cflag(i,j,k-kk-1).isCovered()) {
326 else if (rho_w_avg<0.)
328 for (
int kk=0; kk<ncells_v_vert; ++kk) {
329 if (v_cflag(i,j,k+kk).isCovered()) {
338 amrex::Real interpz = 0.;
339 if (kcell==ncells_v_vert) {
340 interp_v_vert.InterpolateInZ(i,j,k,0,interpz,rho_w_avg);
344 }
else if (kcell==2) {
348 flx_v_arr[2](i,j,k) = rho_w_avg * interpz;
352 flx_v_arr[2](i,j,k) = 0.;
360 amrex::ParallelFor(bxz_grown[0], bxz_grown[1], bxz_grown[2],
362 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
364 if ( w_afrac_x(i,j,k) > 0. )
366 amrex::Real rho_u_avg = 0.5 * (rho_u(i,j,k) + rho_u(i,j,k-1)) * mf_ux_inv(i,j,0);
372 for (
int ii=0; ii<ncells_w_horz; ++ii) {
373 if (w_cflag(i-ii-1,j,k).isCovered()) {
379 else if (rho_u_avg<0.)
381 for (
int ii=0; ii<ncells_w_horz; ++ii) {
382 if (w_cflag(i+ii,j,k).isCovered()) {
391 amrex::Real interpx = 0.;
392 if (icell==ncells_w_horz) {
393 interp_w_horz.InterpolateInX(i,j,k,0,interpx,rho_u_avg);
397 }
else if (icell==2) {
401 flx_w_arr[0](i,j,k) = rho_u_avg * interpx;
405 flx_w_arr[0](i,j,k) = 0.;
409 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
411 if ( w_afrac_y(i,j,k) > 0. )
413 amrex::Real rho_v_avg = 0.5 * (rho_v(i,j,k) + rho_v(i,j,k-1)) * mf_vy_inv(i,j,0);
419 for (
int jj=0; jj<ncells_w_horz; ++jj) {
420 if (w_cflag(i,j-jj-1,k).isCovered()) {
426 else if (rho_v_avg<0.)
428 for (
int jj=0; jj<ncells_w_horz; ++jj) {
429 if (w_cflag(i,j+jj,k).isCovered()) {
438 amrex::Real interpy = 0.;
439 if (jcell==ncells_w_horz) {
440 interp_w_horz.InterpolateInY(i,j,k,0,interpy,rho_v_avg);
444 }
else if (jcell==2) {
448 flx_w_arr[1](i,j,k) = rho_v_avg * interpy;
452 flx_w_arr[1](i,j,k) = 0.;
456 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
458 if ( w_afrac_z(i,j,k) > 0. )
460 if (k==lo_z_face || k==hi_z_face+1) {
462 flx_w_arr[2](i,j,k) = rho_w(i,j,k) * w(i,j,k);
465 amrex::Real rho_w_avg = 0.5 * (rho_w(i,j,k) + rho_w(i,j,k-1));
471 for (
int kk=0; kk<ncells_w_vert; ++kk) {
472 if (w_cflag(i,j,k-kk-1).isCovered()) {
478 else if (rho_w_avg<0.)
480 for (
int kk=0; kk<ncells_w_vert; ++kk) {
481 if (w_cflag(i,j,k+kk).isCovered()) {
490 if (k == hi_z_face || k == lo_z_face+1) {
491 kcell=std::min(kcell,1);
492 }
else if (k == hi_z_face-1 || k == lo_z_face+2 ) {
493 kcell=std::min(kcell,2);
498 amrex::Real interpz = 0.;
499 if (kcell==ncells_w_vert) {
500 interp_w_vert.InterpolateInZ(i,j,k,0,interpz,rho_w_avg);
504 }
else if (kcell==2) {
508 flx_w_arr[2](i,j,k) = rho_w_avg * interpz;
513 flx_w_arr[2](i,j,k) = 0.;
522 template<
typename InterpType_H>
525 const amrex::Vector<amrex::Box>& bxy_grown,
526 const amrex::Vector<amrex::Box>& bxz_grown,
527 const amrex::Array4<const amrex::Real>& rho_u,
528 const amrex::Array4<const amrex::Real>& rho_v,
529 const amrex::Array4<const amrex::Real>& rho_w,
530 const amrex::Array4<const amrex::Real>& u,
531 const amrex::Array4<const amrex::Real>& v,
532 const amrex::Array4<const amrex::Real>& w,
533 const amrex::Array4<const amrex::EBCellFlag>& u_cflag,
534 const amrex::Array4<const amrex::Real>& u_afrac_x,
535 const amrex::Array4<const amrex::Real>& u_afrac_y,
536 const amrex::Array4<const amrex::Real>& u_afrac_z,
537 const amrex::Array4<const amrex::EBCellFlag>& v_cflag,
538 const amrex::Array4<const amrex::Real>& v_afrac_x,
539 const amrex::Array4<const amrex::Real>& v_afrac_y,
540 const amrex::Array4<const amrex::Real>& v_afrac_z,
541 const amrex::Array4<const amrex::EBCellFlag>& w_cflag,
542 const amrex::Array4<const amrex::Real>& w_afrac_x,
543 const amrex::Array4<const amrex::Real>& w_afrac_y,
544 const amrex::Array4<const amrex::Real>& w_afrac_z,
545 const amrex::Array4<const amrex::Real>& mf_ux_inv,
546 const amrex::Array4<const amrex::Real>& mf_vx_inv,
547 const amrex::Array4<const amrex::Real>& mf_uy_inv,
548 const amrex::Array4<const amrex::Real>& mf_vy_inv,
549 const amrex::Real upw_frac_h,
550 const amrex::Real upw_frac_v,
552 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_u_arr,
553 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_v_arr,
554 amrex::GpuArray<amrex::Array4<amrex::Real>, AMREX_SPACEDIM>& flx_w_arr,
555 const int lo_z_face,
const int hi_z_face)
558 EBAdvectionSrcForMomWrapper<InterpType_H,CENTERED2,UPWINDALL>(
559 bxx_grown, bxy_grown, bxz_grown,
560 rho_u, rho_v, rho_w, u, v, w,
561 u_cflag, u_afrac_x, u_afrac_y, u_afrac_z,
562 v_cflag, v_afrac_x, v_afrac_y, v_afrac_z,
563 w_cflag, w_afrac_x, w_afrac_y, w_afrac_z,
564 mf_ux_inv, mf_vx_inv,
565 mf_uy_inv, mf_vy_inv,
566 upw_frac_h, upw_frac_v, vert_adv_type,
567 flx_u_arr, flx_v_arr, flx_w_arr,
568 lo_z_face, hi_z_face);
570 EBAdvectionSrcForMomWrapper<InterpType_H,UPWIND3,UPWINDALL>(
571 bxx_grown, bxy_grown, bxz_grown,
572 rho_u, rho_v, rho_w, u, v, w,
573 u_cflag, u_afrac_x, u_afrac_y, u_afrac_z,
574 v_cflag, v_afrac_x, v_afrac_y, v_afrac_z,
575 w_cflag, w_afrac_x, w_afrac_y, w_afrac_z,
576 mf_ux_inv, mf_vx_inv,
577 mf_uy_inv, mf_vy_inv,
578 upw_frac_h, upw_frac_v, vert_adv_type,
579 flx_u_arr, flx_v_arr, flx_w_arr,
580 lo_z_face, hi_z_face);
582 EBAdvectionSrcForMomWrapper<InterpType_H,CENTERED4,UPWINDALL>(
583 bxx_grown, bxy_grown, bxz_grown,
584 rho_u, rho_v, rho_w, u, v, w,
585 u_cflag, u_afrac_x, u_afrac_y, u_afrac_z,
586 v_cflag, v_afrac_x, v_afrac_y, v_afrac_z,
587 w_cflag, w_afrac_x, w_afrac_y, w_afrac_z,
588 mf_ux_inv, mf_vx_inv,
589 mf_uy_inv, mf_vy_inv,
590 upw_frac_h, upw_frac_v, vert_adv_type,
591 flx_u_arr, flx_v_arr, flx_w_arr,
592 lo_z_face, hi_z_face);
594 EBAdvectionSrcForMomWrapper<InterpType_H,UPWIND5,UPWINDALL>(
595 bxx_grown, bxy_grown, bxz_grown,
596 rho_u, rho_v, rho_w, u, v, w,
597 u_cflag, u_afrac_x, u_afrac_y, u_afrac_z,
598 v_cflag, v_afrac_x, v_afrac_y, v_afrac_z,
599 w_cflag, w_afrac_x, w_afrac_y, w_afrac_z,
600 mf_ux_inv, mf_vx_inv,
601 mf_uy_inv, mf_vy_inv,
602 upw_frac_h, upw_frac_v, vert_adv_type,
603 flx_u_arr, flx_v_arr, flx_w_arr,
604 lo_z_face, hi_z_face);
606 EBAdvectionSrcForMomWrapper<InterpType_H,CENTERED6,UPWINDALL>(
607 bxx_grown, bxy_grown, bxz_grown,
608 rho_u, rho_v, rho_w, u, v, w,
609 u_cflag, u_afrac_x, u_afrac_y, u_afrac_z,
610 v_cflag, v_afrac_x, v_afrac_y, v_afrac_z,
611 w_cflag, w_afrac_x, w_afrac_y, w_afrac_z,
612 mf_ux_inv, mf_vx_inv,
613 mf_uy_inv, mf_vy_inv,
614 upw_frac_h, upw_frac_v, vert_adv_type,
615 flx_u_arr, flx_v_arr, flx_w_arr,
616 lo_z_face, hi_z_face);
618 AMREX_ASSERT_WITH_MESSAGE(
false,
"Unknown advection scheme!");
AdvType
Definition: ERF_IndexDefines.H:221
Definition: ERF_Interpolation_UPW.H:10
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const
Definition: ERF_Interpolation_UPW.H:36
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const
Definition: ERF_Interpolation_UPW.H:54
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const
Definition: ERF_Interpolation_UPW.H:18
Definition: ERF_Interpolation_UPW.H:93
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_UPW.H:155
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_UPW.H:129
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const
Definition: ERF_Interpolation_UPW.H:103