Function for computing the scalar RHS for diffusion operator with terrain-fitted coordinates.
88 const Real explicit_fac = 1.0 - 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;
115 rhoAlpha += 0.5 * ( mu_turb(i , j, k,
d_eddy_diff_idx[prim_scal_index])
124 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
126 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
127 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
128 Real GradCz = 0.5 * ( 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;
147 rhoAlpha += 0.5 * ( mu_turb(i, j , k,
d_eddy_diff_idy[prim_scal_index])
155 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
157 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
158 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
159 Real GradCz = 0.5 * ( 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;
178 rhoAlpha += 0.5 * ( mu_turb(i, j, k ,
d_eddy_diff_idz[prim_scal_index])
191 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
193 if (ext_dir_on_zlo) {
198 Real idz0 = 1.0 / dz0;
199 Real f = (dz1 / dz0) + 2.0;
201 Real c3 = 2.0 / (f - f2);
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) {
212 Real idz0 = 1.0 / dz0;
213 Real f = (dz1 / dz0) + 2.0;
215 Real c3 = 2.0 / (f - f2);
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) );
227 zflux(i,j,k) = hfx_z(i,j,0);
228 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
229 zflux(i,j,k) = qfx1_z(i,j,0);
231 zflux(i,j,k) = -rhoAlpha * GradCz;
235 if (!SurfLayer_on_zlo) {
236 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
239 if (!SurfLayer_on_zlo) {
240 qfx1_z(i,j,k) = zflux(i,j,k);
243 qfx2_z(i,j,k) = zflux(i,j,k);
248 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
250 const int prim_index = qty_index - 1;
261 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
263 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
264 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
265 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
266 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
267 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
270 xflux(i,j,k) = hfx_x(i,j,0);
271 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
272 xflux(i,j,k) = qfx1_x(i,j,0);
274 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
277 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
279 const int prim_index = qty_index - 1;
290 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
292 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
293 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
294 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
295 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
296 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
299 yflux(i,j,k) = hfx_y(i,j,0);
300 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
301 yflux(i,j,k) = qfx1_y(i,j,0);
303 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
306 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
308 const int prim_index = qty_index - 1;
324 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
326 if (ext_dir_on_zlo) {
331 Real idz0 = 1.0 / dz0;
332 Real f = (dz1 / dz0) + 2.0;
334 Real c3 = 2.0 / (f - f2);
337 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
338 +
c2 * cell_prim(i, j, k , prim_index)
339 + c3 * cell_prim(i, j, k+1, prim_index) );
340 }
else if (ext_dir_on_zhi) {
345 Real idz0 = 1.0 / dz0;
346 Real f = (dz1 / dz0) + 2.0;
348 Real c3 = 2.0 / (f - f2);
351 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
352 +
c2 * cell_prim(i, j, k-1, prim_index)
353 + c3 * cell_prim(i, j, k-2, prim_index) ) );
356 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
360 zflux(i,j,k) = hfx_z(i,j,0);
361 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
362 zflux(i,j,k) = qfx1_z(i,j,0);
364 zflux(i,j,k) = -rhoAlpha * GradCz;
368 if (!SurfLayer_on_zlo) {
369 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
372 if (!SurfLayer_on_zlo) {
373 qfx1_z(i,j,k) = zflux(i,j,k);
376 qfx2_z(i,j,k) = zflux(i,j,k);
382 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
384 const int prim_index = qty_index - 1;
394 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
396 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
397 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
398 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
399 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
400 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
403 xflux(i,j,k) = hfx_x(i,j,0);
404 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
405 xflux(i,j,k) = qfx1_x(i,j,0);
407 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
410 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
412 const int prim_index = qty_index - 1;
422 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
424 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
425 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
426 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
427 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
428 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
431 yflux(i,j,k) = hfx_y(i,j,0);
432 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
433 yflux(i,j,k) = qfx1_y(i,j,0);
435 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
438 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
440 const int prim_index = qty_index - 1;
455 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
457 if (ext_dir_on_zlo) {
462 Real idz0 = 1.0 / dz0;
463 Real f = (dz1 / dz0) + 2.0;
465 Real c3 = 2.0 / (f - f2);
468 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
469 +
c2 * cell_prim(i, j, k , prim_index)
470 + c3 * cell_prim(i, j, k+1, prim_index) );
471 }
else if (ext_dir_on_zhi) {
476 Real idz0 = 1.0 / dz0;
477 Real f = (dz1 / dz0) + 2.0;
479 Real c3 = 2.0 / (f - f2);
482 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
483 +
c2 * cell_prim(i, j, k-1, prim_index)
484 + c3 * cell_prim(i, j, k-2, prim_index) ) );
487 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
491 zflux(i,j,k) = hfx_z(i,j,0);
492 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
493 zflux(i,j,k) = qfx1_z(i,j,0);
495 zflux(i,j,k) = -rhoAlpha * GradCz;
499 if (!SurfLayer_on_zlo) {
500 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
503 if (!SurfLayer_on_zlo) {
504 qfx1_z(i,j,k) = zflux(i,j,k);
507 qfx2_z(i,j,k) = zflux(i,j,k);
512 ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
514 const int prim_index = qty_index - 1;
523 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
525 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
526 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
527 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
528 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
529 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
532 xflux(i,j,k) = hfx_x(i,j,0);
533 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
534 xflux(i,j,k) = qfx1_x(i,j,0);
536 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
539 ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
541 const int prim_index = qty_index - 1;
550 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
552 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
553 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
554 Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
555 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
556 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
559 yflux(i,j,k) = hfx_y(i,j,0);
560 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
561 yflux(i,j,k) = qfx1_y(i,j,0);
563 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
566 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
568 const int prim_index = qty_index - 1;
583 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
585 if (ext_dir_on_zlo) {
590 Real idz0 = 1.0 / dz0;
591 Real f = (dz1 / dz0) + 2.0;
593 Real c3 = 2.0 / (f - f2);
596 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
597 +
c2 * cell_prim(i, j, k , prim_index)
598 + c3 * cell_prim(i, j, k+1, prim_index) );
599 }
else if (ext_dir_on_zhi) {
604 Real idz0 = 1.0 / dz0;
605 Real f = (dz1 / dz0) + 2.0;
607 Real c3 = 2.0 / (f - f2);
610 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
611 +
c2 * cell_prim(i, j, k-1, prim_index)
612 + c3 * cell_prim(i, j, k-2, prim_index) ) );
615 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
619 zflux(i,j,k) = hfx_z(i,j,0);
620 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
621 zflux(i,j,k) = qfx1_z(i,j,0);
623 zflux(i,j,k) = -rhoAlpha * GradCz;
627 if (!SurfLayer_on_zlo) {
628 hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
631 if (!SurfLayer_on_zlo) {
632 qfx1_z(i,j,k) = zflux(i,j,k);
635 qfx2_z(i,j,k) = zflux(i,j,k);
647 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
655 Real xfluxbar_lo, yfluxbar_lo;
657 Real xfluxlo = 0.5 * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
658 Real xfluxhi = 0.5 * ( xflux(i,j,k+1) + xflux(i+1,j,k+1) );
659 xfluxbar_lo = 1.5*xfluxlo - 0.5*xfluxhi;
661 Real yfluxlo = 0.5 * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
662 Real yfluxhi = 0.5 * ( yflux(i,j,k+1) + yflux(i,j+1,k+1) );
663 yfluxbar_lo = 1.5*yfluxlo - 0.5*yfluxhi;
665 xfluxbar_lo = 0.25 * ( xflux(i,j,k ) + xflux(i+1,j ,k )
666 + xflux(i,j,k-1) + xflux(i+1,j ,k-1) );
667 yfluxbar_lo = 0.25 * ( yflux(i,j,k ) + yflux(i ,j+1,k )
668 + yflux(i,j,k-1) + yflux(i ,j+1,k-1) );
671 Real xfluxbar_hi, yfluxbar_hi;
673 Real xfluxlo = 0.5 * ( xflux(i,j,k-1) + xflux(i+1,j,k-1) );
674 Real xfluxhi = 0.5 * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
675 xfluxbar_hi = 1.5*xfluxhi - 0.5*xfluxlo;
677 Real yfluxlo = 0.5 * ( yflux(i,j,k-1) + yflux(i,j+1,k-1) );
678 Real yfluxhi = 0.5 * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
679 yfluxbar_hi = 1.5*yfluxhi - 0.5*yfluxlo;
681 xfluxbar_hi = 0.25 * ( xflux(i,j,k+1) + xflux(i+1,j ,k+1)
682 + xflux(i,j,k ) + xflux(i+1,j ,k ) );
683 yfluxbar_hi = 0.25 * ( yflux(i,j,k+1) + yflux(i ,j+1,k+1)
684 + yflux(i,j,k ) + yflux(i ,j+1,k ) );
688 Real zflux_lo = explicit_fac * zflux(i,j,k )
689 - met_h_xi_lo * mf_mx(i,j,0) * xfluxbar_lo
690 - met_h_eta_lo * mf_my(i,j,0) * yfluxbar_lo;
691 Real zflux_hi = explicit_fac * zflux(i,j,k+1)
692 - met_h_xi_hi * mf_mx(i,j,0) * xfluxbar_hi
693 - met_h_eta_hi * mf_my(i,j,0) * yfluxbar_hi;
695 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
696 Real stateContrib = ( xflux(i+1,j ,k ) * ax(i+1,j,k) / mf_uy(i+1,j,0)
697 -xflux(i ,j ,k ) * ax(i ,j,k) / mf_uy(i ,j,0) ) *
dx_inv * mfsq
698 +( yflux(i ,j+1,k ) * ay(i,j+1,k) / mf_vx(i,j+1,0)
699 -yflux(i ,j ,k ) * ay(i,j ,k) / mf_vx(i,j ,0) ) *
dy_inv * mfsq
700 +( zflux_hi - zflux_lo) * dz_inv;
702 stateContrib /= detJ(i,j,k);
704 cell_rhs(i,j,k,qty_index) -= stateContrib;
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
const int bc_comp
Definition: ERF_Implicit.H:8
#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:52
int * d_eddy_diff_idx
Definition: ERF_SetupDiff.H:47
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:48
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:8
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
bool l_consA
Definition: ERF_SetupVertDiff.H:7
int * d_eddy_diff_idz
Definition: ERF_SetupVertDiff.H:107
Real * d_alpha_eff
Definition: ERF_SetupVertDiff.H:104
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:115
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:182
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:374
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:196
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 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:209
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:80
@ ext_dir
Definition: ERF_IndexDefines.H:209
@ ext_dir_prim
Definition: ERF_IndexDefines.H:211
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