Function for computing the scalar RHS for diffusion operator with terrain-fitted coordinates.
88 const Real explicit_fac =
one - implicit_fac;
91 Real l_abs_g = std::abs(grav_gpu[2]);
93 const Real dz_inv = cellSizeInv[2];
97 Box xbx_g1(
xbx); Box ybx_g1(
ybx);
98 if (xbx_g1.smallEnd(2) !=
dom_lo.z) xbx_g1.growLo(2,1);
99 if (ybx_g1.smallEnd(2) !=
dom_lo.z) ybx_g1.growLo(2,1);
100 if (xbx_g1.bigEnd(2) !=
dom_hi.z) xbx_g1.growHi(2,1);
101 if (ybx_g1.bigEnd(2) !=
dom_hi.z) ybx_g1.growHi(2,1);
103 for (
int n(0); n<num_comp; ++n) {
104 const int qty_index = start_comp + n;
108 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
110 const int prim_index = qty_index - 1;
124 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
126 Real idz_hi =
one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
127 Real idz_lo =
one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
128 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
129 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
130 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
133 xflux(i,j,k) = hfx_x(i,j,0);
134 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
135 xflux(i,j,k) = qfx1_x(i,j,0);
137 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
140 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
142 const int prim_index = qty_index - 1;
155 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
157 Real idz_hi =
one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
158 Real idz_lo =
one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
159 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
160 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
161 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
164 yflux(i,j,k) = hfx_y(i,j,0);
165 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
166 yflux(i,j,k) = qfx1_y(i,j,0);
168 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
171 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
173 const int prim_index = qty_index - 1;
191 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
193 if (ext_dir_on_zlo) {
204 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
205 +
c2 * cell_prim(i, j, k , prim_index)
206 + c3 * cell_prim(i, j, k+1, prim_index) );
207 }
else if (ext_dir_on_zhi) {
218 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
219 +
c2 * cell_prim(i, j, k-1, prim_index)
220 + c3 * cell_prim(i, j, k-2, prim_index) ) );
223 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
226 if (SurfLayer_on_zlo) {
228 zflux(i,j,k) = hfx_z(i,j,0);
230 zflux(i,j,k) = qfx1_z(i,j,0);
235 zflux(i,j,k) = -rhoAlpha * GradCz;
239 if (!SurfLayer_on_zlo) {
240 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
243 if (!SurfLayer_on_zlo) {
244 qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
247 qfx2_z(i,j,k) = zflux(i,j,k);
252 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
254 const int prim_index = qty_index - 1;
265 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
267 Real idz_hi =
one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
268 Real idz_lo =
one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
269 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
270 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
271 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
274 xflux(i,j,k) = hfx_x(i,j,0);
275 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
276 xflux(i,j,k) = qfx1_x(i,j,0);
278 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
281 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
283 const int prim_index = qty_index - 1;
294 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
296 Real idz_hi =
one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
297 Real idz_lo =
one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
298 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
299 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
300 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
303 yflux(i,j,k) = hfx_y(i,j,0);
304 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
305 yflux(i,j,k) = qfx1_y(i,j,0);
307 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
310 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
312 const int prim_index = qty_index - 1;
328 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
330 if (ext_dir_on_zlo) {
341 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
342 +
c2 * cell_prim(i, j, k , prim_index)
343 + c3 * cell_prim(i, j, k+1, prim_index) );
344 }
else if (ext_dir_on_zhi) {
355 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
356 +
c2 * cell_prim(i, j, k-1, prim_index)
357 + c3 * cell_prim(i, j, k-2, prim_index) ) );
360 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
363 if (SurfLayer_on_zlo) {
365 zflux(i,j,k) = hfx_z(i,j,0);
367 zflux(i,j,k) = qfx1_z(i,j,0);
372 zflux(i,j,k) = -rhoAlpha * GradCz;
376 if (!SurfLayer_on_zlo) {
377 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
380 if (!SurfLayer_on_zlo) {
381 qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
384 qfx2_z(i,j,k) = zflux(i,j,k);
390 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
392 const int prim_index = qty_index - 1;
402 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
404 Real idz_hi =
one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
405 Real idz_lo =
one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
406 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
407 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
408 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
411 xflux(i,j,k) = hfx_x(i,j,0);
412 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
413 xflux(i,j,k) = qfx1_x(i,j,0);
415 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
418 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
420 const int prim_index = qty_index - 1;
430 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
432 Real idz_hi =
one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
433 Real idz_lo =
one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
434 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
435 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
436 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
439 yflux(i,j,k) = hfx_y(i,j,0);
440 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
441 yflux(i,j,k) = qfx1_y(i,j,0);
443 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
446 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
448 const int prim_index = qty_index - 1;
463 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
465 if (ext_dir_on_zlo) {
476 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
477 +
c2 * cell_prim(i, j, k , prim_index)
478 + c3 * cell_prim(i, j, k+1, prim_index) );
479 }
else if (ext_dir_on_zhi) {
490 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
491 +
c2 * cell_prim(i, j, k-1, prim_index)
492 + c3 * cell_prim(i, j, k-2, prim_index) ) );
495 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
498 if (SurfLayer_on_zlo) {
500 zflux(i,j,k) = hfx_z(i,j,0);
502 zflux(i,j,k) = qfx1_z(i,j,0);
507 zflux(i,j,k) = -rhoAlpha * GradCz;
511 if (!SurfLayer_on_zlo) {
512 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
515 if (!SurfLayer_on_zlo) {
516 qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
519 qfx2_z(i,j,k) = zflux(i,j,k);
524 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
526 const int prim_index = qty_index - 1;
535 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
537 Real idz_hi =
one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
538 Real idz_lo =
one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
539 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
540 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
541 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
544 xflux(i,j,k) = hfx_x(i,j,0);
545 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
546 xflux(i,j,k) = qfx1_x(i,j,0);
548 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
551 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
553 const int prim_index = qty_index - 1;
562 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
564 Real idz_hi =
one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
565 Real idz_lo =
one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
566 Real GradCz =
myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
567 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
568 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
571 yflux(i,j,k) = hfx_y(i,j,0);
572 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
573 yflux(i,j,k) = qfx1_y(i,j,0);
575 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
578 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
580 const int prim_index = qty_index - 1;
595 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
597 if (ext_dir_on_zlo) {
608 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
609 +
c2 * cell_prim(i, j, k , prim_index)
610 + c3 * cell_prim(i, j, k+1, prim_index) );
611 }
else if (ext_dir_on_zhi) {
622 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
623 +
c2 * cell_prim(i, j, k-1, prim_index)
624 + c3 * cell_prim(i, j, k-2, prim_index) ) );
627 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
630 if (SurfLayer_on_zlo) {
632 zflux(i,j,k) = hfx_z(i,j,0);
634 zflux(i,j,k) = qfx1_z(i,j,0);
639 zflux(i,j,k) = -rhoAlpha * GradCz;
643 if (!SurfLayer_on_zlo) {
644 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
647 if (!SurfLayer_on_zlo) {
648 qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
651 qfx2_z(i,j,k) = zflux(i,j,k);
660 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
662 zflux(i,j,k) *= explicit_fac;
673 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
675 Real xfluxbar_lo, yfluxbar_lo;
679 Real xfluxlo =
myhalf * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
680 Real xfluxhi =
myhalf * ( xflux(i,j,k+1) + xflux(i+1,j,k+1) );
681 xfluxbar_lo =
Real(1.5)*xfluxlo -
myhalf*xfluxhi;
683 Real yfluxlo =
myhalf * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
684 Real yfluxhi =
myhalf * ( yflux(i,j,k+1) + yflux(i,j+1,k+1) );
685 yfluxbar_lo =
Real(1.5)*yfluxlo -
myhalf*yfluxhi;
687 xfluxbar_lo =
fourth * ( xflux(i,j,k ) + xflux(i+1,j ,k )
688 + xflux(i,j,k-1) + xflux(i+1,j ,k-1) );
689 yfluxbar_lo =
fourth * ( yflux(i,j,k ) + yflux(i ,j+1,k )
690 + yflux(i,j,k-1) + yflux(i ,j+1,k-1) );
693 Real xfluxbar_hi, yfluxbar_hi;
697 Real xfluxlo =
myhalf * ( xflux(i,j,k-1) + xflux(i+1,j,k-1) );
698 Real xfluxhi =
myhalf * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
699 xfluxbar_hi =
Real(1.5)*xfluxhi -
myhalf*xfluxlo;
701 Real yfluxlo =
myhalf * ( yflux(i,j,k-1) + yflux(i,j+1,k-1) );
702 Real yfluxhi =
myhalf * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
703 yfluxbar_hi =
Real(1.5)*yfluxhi -
myhalf*yfluxlo;
705 xfluxbar_hi =
fourth * ( xflux(i,j,k+1) + xflux(i+1,j ,k+1)
706 + xflux(i,j,k ) + xflux(i+1,j ,k ) );
707 yfluxbar_hi =
fourth * ( yflux(i,j,k+1) + yflux(i ,j+1,k+1)
708 + yflux(i,j,k ) + yflux(i ,j+1,k ) );
713 if ( use_SurfLayer &&
719 zflux_lo = zflux(i,j,k )
720 - met_h_xi_lo * mf_mx(i,j,0) * xfluxbar_lo
721 - met_h_eta_lo * mf_my(i,j,0) * yfluxbar_lo;
723 Real zflux_hi = zflux(i,j,k+1)
724 - met_h_xi_hi * mf_mx(i,j,0) * xfluxbar_hi
725 - met_h_eta_hi * mf_my(i,j,0) * yfluxbar_hi;
727 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
728 Real stateContrib = ( xflux(i+1,j ,k ) * ax(i+1,j,k) / mf_uy(i+1,j,0)
729 -xflux(i ,j ,k ) * ax(i ,j,k) / mf_uy(i ,j,0) ) *
dx_inv * mfsq
730 +( yflux(i ,j+1,k ) * ay(i,j+1,k) / mf_vx(i,j+1,0)
731 -yflux(i ,j ,k ) * ay(i,j ,k) / mf_vx(i,j ,0) ) *
dy_inv * mfsq
732 +( zflux_hi - zflux_lo) * dz_inv;
734 stateContrib /= detJ(i,j,k);
736 cell_rhs(i,j,k,qty_index) -= stateContrib;
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
void DiffusionSrcForState_T(const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &rotate, const Array4< const Real > &u, const Array4< const Real > &v, const Array4< const Real > &cell_data, const Array4< const Real > &cell_prim, const Array4< Real > &cell_rhs, const Array4< Real > &xflux, const Array4< Real > &yflux, const Array4< Real > &zflux, const Array4< const Real > &z_nd, const Array4< const Real > &z_cc, const Array4< const Real > &ax, const Array4< const Real > &ay, const Array4< const Real > &, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_ux, const Array4< const Real > &mf_vx, const Array4< const Real > &mf_my, const Array4< const Real > &mf_uy, const Array4< const Real > &mf_vy, Array4< Real > &hfx_x, Array4< Real > &hfx_y, Array4< Real > &hfx_z, Array4< Real > &qfx1_x, Array4< Real > &qfx1_y, Array4< Real > &qfx1_z, Array4< Real > &qfx2_z, Array4< Real > &diss, const Array4< const Real > &mu_turb, const SolverChoice &solverChoice, const int level, const Array4< const Real > &tm_arr, const GpuArray< Real, AMREX_SPACEDIM > grav_gpu, const BCRec *bc_ptr, const bool use_SurfLayer, const Real implicit_fac)
Definition: ERF_DiffusionSrcForState_T.cpp:44
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define PrimScalar_comp
Definition: ERF_IndexDefines.H:57
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);})
int * d_eddy_diff_idx
Definition: ERF_SetupDiff.H:45
const Box zbx
Definition: ERF_SetupDiff.H:9
const Real dx_inv
Definition: ERF_SetupDiff.H:4
const Real dy_inv
Definition: ERF_SetupDiff.H:5
int * d_eddy_diff_idy
Definition: ERF_SetupDiff.H:46
const Box xbx
Definition: ERF_SetupDiff.H:7
const Box ybx
Definition: ERF_SetupDiff.H:8
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
bool l_turb
Definition: ERF_SetupVertDiff.H:9
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
bool l_consA
Definition: ERF_SetupVertDiff.H:8
int * d_eddy_diff_idz
Definition: ERF_SetupVertDiff.H:99
Real * d_alpha_eff
Definition: ERF_SetupVertDiff.H:98
amrex::Real Real
Definition: ERF_ShocInterface.H:19
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:117
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface(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:184
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:376
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtKface(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:198
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:170
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtKface(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:211
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:90
@ ext_dir
Definition: ERF_IndexDefines.H:243
@ ext_dir_prim
Definition: ERF_IndexDefines.H:246
real(c_double), parameter c2
Definition: ERF_module_model_constants.F90:35
real(c_double), private c1
Definition: ERF_module_mp_morr_two_moment.F90:212