Function for computing the scalar RHS for diffusion operator with terrain-fitted coordinates.
88 const Real dz_inv = cellSizeInv[2];
92 for (
int n(0); n<num_comp; ++n) {
93 const int qty_index = start_comp + n;
97 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
99 const int prim_index = qty_index - 1;
102 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
103 Real rhoAlpha = rhoFace *
d_alpha_eff[prim_scal_index];
104 rhoAlpha += 0.5 * ( mu_turb(i , j, k,
d_eddy_diff_idx[prim_scal_index])
113 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
115 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
116 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
117 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
118 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
119 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
122 xflux(i,j,k) = hfx_x(i,j,0);
123 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
124 xflux(i,j,k) = qfx1_x(i,j,0);
126 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
129 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
131 const int prim_index = qty_index - 1;
134 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
135 Real rhoAlpha = rhoFace *
d_alpha_eff[prim_scal_index];
136 rhoAlpha += 0.5 * ( mu_turb(i, j , k,
d_eddy_diff_idy[prim_scal_index])
144 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
146 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
147 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
148 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
149 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
150 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
153 yflux(i,j,k) = hfx_y(i,j,0);
154 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
155 yflux(i,j,k) = qfx1_y(i,j,0);
157 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
160 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
162 const int prim_index = qty_index - 1;
165 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
166 Real rhoAlpha = rhoFace *
d_alpha_eff[prim_scal_index];
167 rhoAlpha += 0.5 * ( mu_turb(i, j, k ,
d_eddy_diff_idz[prim_scal_index])
182 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
184 if (ext_dir_on_zlo) {
189 Real idz0 = 1.0 / dz0;
190 Real f = (dz1 / dz0) + 2.0;
192 Real c3 = 2.0 / (f - f2);
194 Real
c1 = -(1.0-f2)*c3;
195 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
196 +
c2 * cell_prim(i, j, k , prim_index)
197 + c3 * cell_prim(i, j, k+1, prim_index) );
198 }
else if (ext_dir_on_zhi) {
203 Real idz0 = 1.0 / dz0;
204 Real f = (dz1 / dz0) + 2.0;
206 Real c3 = 2.0 / (f - f2);
208 Real
c1 = -(1.0-f2)*c3;
209 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
210 +
c2 * cell_prim(i, j, k-1, prim_index)
211 + c3 * cell_prim(i, j, k-2, prim_index) ) );
213 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
217 zflux(i,j,k) = hfx_z(i,j,0);
218 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
219 zflux(i,j,k) = qfx1_z(i,j,0);
221 zflux(i,j,k) = -rhoAlpha * GradCz;
225 if (!SurfLayer_on_zlo) {
226 hfx_z(i,j,k) = zflux(i,j,k);
229 if (!SurfLayer_on_zlo) {
230 qfx1_z(i,j,k) = zflux(i,j,k);
233 qfx2_z(i,j,k) = zflux(i,j,k);
238 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
240 const int prim_index = qty_index - 1;
251 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
253 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
254 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
255 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
256 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
257 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
260 xflux(i,j,k) = hfx_x(i,j,0);
261 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
262 xflux(i,j,k) = qfx1_x(i,j,0);
264 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
267 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
269 const int prim_index = qty_index - 1;
280 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
282 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
283 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
284 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
285 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
286 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
289 yflux(i,j,k) = hfx_y(i,j,0);
290 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
291 yflux(i,j,k) = qfx1_y(i,j,0);
293 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
296 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
298 const int prim_index = qty_index - 1;
317 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
319 if (ext_dir_on_zlo) {
324 Real idz0 = 1.0 / dz0;
325 Real f = (dz1 / dz0) + 2.0;
327 Real c3 = 2.0 / (f - f2);
329 Real
c1 = -(1.0-f2)*c3;
330 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
331 +
c2 * cell_prim(i, j, k , prim_index)
332 + c3 * cell_prim(i, j, k+1, prim_index) );
333 }
else if (ext_dir_on_zhi) {
338 Real idz0 = 1.0 / dz0;
339 Real f = (dz1 / dz0) + 2.0;
341 Real c3 = 2.0 / (f - f2);
343 Real
c1 = -(1.0-f2)*c3;
344 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
345 +
c2 * cell_prim(i, j, k-1, prim_index)
346 + c3 * cell_prim(i, j, k-2, prim_index) ) );
348 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
352 zflux(i,j,k) = hfx_z(i,j,0);
353 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
354 zflux(i,j,k) = qfx1_z(i,j,0);
356 zflux(i,j,k) = -rhoAlpha * GradCz;
360 if (!SurfLayer_on_zlo) {
361 hfx_z(i,j,k) = zflux(i,j,k);
364 if (!SurfLayer_on_zlo) {
365 qfx1_z(i,j,k) = zflux(i,j,k);
368 qfx2_z(i,j,k) = zflux(i,j,k);
373 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
375 const int prim_index = qty_index - 1;
377 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
385 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
387 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
388 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
389 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
390 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
391 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
394 xflux(i,j,k) = hfx_x(i,j,0);
395 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
396 xflux(i,j,k) = qfx1_x(i,j,0);
398 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
401 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
403 const int prim_index = qty_index - 1;
405 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
413 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
415 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
416 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
417 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
418 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
419 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
422 yflux(i,j,k) = hfx_y(i,j,0);
423 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
424 yflux(i,j,k) = qfx1_y(i,j,0);
426 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
429 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
431 const int prim_index = qty_index - 1;
433 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
448 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
450 if (ext_dir_on_zlo) {
455 Real idz0 = 1.0 / dz0;
456 Real f = (dz1 / dz0) + 2.0;
458 Real c3 = 2.0 / (f - f2);
460 Real
c1 = -(1.0-f2)*c3;
461 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
462 +
c2 * cell_prim(i, j, k , prim_index)
463 + c3 * cell_prim(i, j, k+1, prim_index) );
464 }
else if (ext_dir_on_zhi) {
469 Real idz0 = 1.0 / dz0;
470 Real f = (dz1 / dz0) + 2.0;
472 Real c3 = 2.0 / (f - f2);
474 Real
c1 = -(1.0-f2)*c3;
475 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
476 +
c2 * cell_prim(i, j, k-1, prim_index)
477 + c3 * cell_prim(i, j, k-2, prim_index) ) );
479 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
483 zflux(i,j,k) = hfx_z(i,j,0);
484 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
485 zflux(i,j,k) = qfx1_z(i,j,0);
487 zflux(i,j,k) = -rhoAlpha * GradCz;
491 if (!SurfLayer_on_zlo) {
492 hfx_z(i,j,k) = zflux(i,j,k);
495 if (!SurfLayer_on_zlo) {
496 qfx1_z(i,j,k) = zflux(i,j,k);
499 qfx2_z(i,j,k) = zflux(i,j,k);
504 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
506 const int prim_index = qty_index - 1;
515 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
517 Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
518 Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
519 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
520 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
521 Real GradCx =
dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
524 xflux(i,j,k) = hfx_x(i,j,0);
525 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
526 xflux(i,j,k) = qfx1_x(i,j,0);
528 xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
531 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
533 const int prim_index = qty_index - 1;
542 bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k ==
dom_lo.z);
544 Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
545 Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
546 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
547 - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
548 Real GradCy =
dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
551 yflux(i,j,k) = hfx_y(i,j,0);
552 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
553 yflux(i,j,k) = qfx1_y(i,j,0);
555 yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
558 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
560 const int prim_index = qty_index - 1;
577 bool SurfLayer_on_zlo = ( use_SurfLayer && k ==
dom_lo.z);
579 if (ext_dir_on_zlo) {
584 Real idz0 = 1.0 / dz0;
585 Real f = (dz1 / dz0) + 2.0;
587 Real c3 = 2.0 / (f - f2);
589 Real
c1 = -(1.0-f2)*c3;
590 GradCz = idz0 * (
c1 * cell_prim(i, j, k-1, prim_index)
591 +
c2 * cell_prim(i, j, k , prim_index)
592 + c3 * cell_prim(i, j, k+1, prim_index) );
593 }
else if (ext_dir_on_zhi) {
598 Real idz0 = 1.0 / dz0;
599 Real f = (dz1 / dz0) + 2.0;
601 Real c3 = 2.0 / (f - f2);
603 Real
c1 = -(1.0-f2)*c3;
604 GradCz = idz0 * ( -(
c1 * cell_prim(i, j, k , prim_index)
605 +
c2 * cell_prim(i, j, k-1, prim_index)
606 + c3 * cell_prim(i, j, k-2, prim_index) ) );
608 GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
612 zflux(i,j,k) = hfx_z(i,j,0);
613 }
else if (SurfLayer_on_zlo && (qty_index ==
RhoQ1_comp)) {
614 zflux(i,j,k) = qfx1_z(i,j,0);
616 zflux(i,j,k) = -rhoAlpha * GradCz;
620 if (!SurfLayer_on_zlo) {
621 hfx_z(i,j,k) = zflux(i,j,k);
624 if (!SurfLayer_on_zlo) {
625 qfx1_z(i,j,k) = zflux(i,j,k);
628 qfx2_z(i,j,k) = zflux(i,j,k);
637 Box planexy =
zbx; planexy.setBig(2, planexy.smallEnd(2) );
638 int k_lo =
zbx.smallEnd(2);
int k_hi =
zbx.bigEnd(2);
639 zbx3.growLo(2,-1); zbx3.growHi(2,-1);
640 ParallelFor(planexy, [=] AMREX_GPU_DEVICE (
int i,
int j,
int ) noexcept
642 Real met_h_xi,met_h_eta;
648 Real xfluxlo = 0.5 * ( xflux(i , j , k_lo ) + xflux(i+1, j , k_lo ) );
649 Real xfluxhi = 0.5 * ( xflux(i , j , k_lo+1) + xflux(i+1, j , k_lo+1) );
650 Real xfluxbar = 1.5*xfluxlo - 0.5*xfluxhi;
652 Real yfluxlo = 0.5 * ( yflux(i , j , k_lo ) + yflux(i , j+1, k_lo ) );
653 Real yfluxhi = 0.5 * ( yflux(i , j , k_lo+1) + yflux(i , j+1, k_lo+1) );
654 Real yfluxbar = 1.5*yfluxlo - 0.5*yfluxhi;
656 zflux(i,j,k_lo) -= met_h_xi*mf_mx(i,j,0)*xfluxbar + met_h_eta*mf_my(i,j,0)*yfluxbar;
663 Real xfluxlo = 0.5 * ( xflux(i , j , k_hi-2) + xflux(i+1, j , k_hi-2) );
664 Real xfluxhi = 0.5 * ( xflux(i , j , k_hi-1) + xflux(i+1, j , k_hi-1) );
665 Real xfluxbar = 1.5*xfluxhi - 0.5*xfluxlo;
667 Real yfluxlo = 0.5 * ( yflux(i , j , k_hi-2) + yflux(i , j+1, k_hi-2) );
668 Real yfluxhi = 0.5 * ( yflux(i , j , k_hi-1) + yflux(i , j+1, k_hi-1) );
669 Real yfluxbar = 1.5*yfluxhi - 0.5*yfluxlo;
671 zflux(i,j,k_hi) -= met_h_xi*mf_mx(i,j,0)*xfluxbar + met_h_eta*mf_my(i,j,0)*yfluxbar;
676 ParallelFor(zbx3, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
678 Real met_h_xi,met_h_eta;
682 Real xfluxbar = 0.25 * ( xflux(i , j , k ) + xflux(i+1, j , k )
683 + xflux(i , j , k-1) + xflux(i+1, j , k-1) );
684 Real yfluxbar = 0.25 * ( yflux(i , j , k ) + yflux(i , j+1, k )
685 + yflux(i , j , k-1) + yflux(i , j+1, k-1) );
687 zflux(i,j,k) -= met_h_xi*mf_mx(i,j,0)*xfluxbar + met_h_eta*mf_my(i,j,0)*yfluxbar;
690 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
692 xflux(i,j,k) *= ax(i,j,k)/mf_uy(i,j,0);
694 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
696 yflux(i,j,k) *= ay(i,j,k)/mf_vx(i,j,0);
702 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
704 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
705 Real stateContrib = (xflux(i+1,j ,k ) - xflux(i, j, k)) *
dx_inv * mfsq
706 +(yflux(i ,j+1,k ) - yflux(i, j, k)) *
dy_inv * mfsq
707 +(zflux(i ,j ,k+1) - zflux(i, j, k)) * dz_inv;
709 stateContrib /= detJ(i,j,k);
711 cell_rhs(i,j,k,qty_index) -= stateContrib;
const auto & dom_hi
Definition: ERF_DiffSetup.H:10
bool l_turb
Definition: ERF_DiffSetup.H:19
const auto & dom_lo
Definition: ERF_DiffSetup.H:9
bool l_consA
Definition: ERF_DiffSetup.H:18
int * d_eddy_diff_idx
Definition: ERF_DiffSetup.H:120
const Box zbx
Definition: ERF_DiffSetup.H:23
int * d_eddy_diff_idz
Definition: ERF_DiffSetup.H:122
const Real dx_inv
Definition: ERF_DiffSetup.H:6
const Real dy_inv
Definition: ERF_DiffSetup.H:7
int * d_eddy_diff_idy
Definition: ERF_DiffSetup.H:121
const Box xbx
Definition: ERF_DiffSetup.H:21
Real * d_alpha_eff
Definition: ERF_DiffSetup.H:119
const Box ybx
Definition: ERF_DiffSetup.H:22
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)
Definition: ERF_DiffusionSrcForState_T.cpp:43
#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
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:110
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:377
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:197
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:167
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: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:211