Function for computing the scalar RHS for diffusion operator without terrain.
73 amrex::ignore_unused(use_most);
75 const Real dx_inv = cellSizeInv[0];
76 const Real dy_inv = cellSizeInv[1];
77 const Real dz_inv = cellSizeInv[2];
79 const auto& dom_lo = lbound(domain);
80 const auto& dom_hi = ubound(domain);
82 Real l_inv_theta0 = 1.0 / turbChoice.
theta_ref;
83 Real l_abs_g = std::abs(grav_gpu[2]);
85 bool l_use_keqn = ( (turbChoice.
les_type == LESType::Deardorff) ||
86 (turbChoice.
rans_type == RANSType::kEqn) );
87 bool l_use_mynn = (turbChoice.
pbl_type == PBLType::MYNN25);
90 bool l_turb = ( (turbChoice.
les_type == LESType::Smagorinsky) ||
91 (turbChoice.
les_type == LESType::Deardorff ) ||
92 (turbChoice.
rans_type == RANSType::kEqn ) ||
93 (turbChoice.
pbl_type == PBLType::MYNN25 ) ||
94 (turbChoice.
pbl_type == PBLType::YSU ) );
96 const Box xbx = surroundingNodes(bx,0);
97 const Box ybx = surroundingNodes(bx,1);
98 const Box zbx = surroundingNodes(bx,2);
100 const int end_comp = start_comp + num_comp - 1;
181 Gpu::AsyncVector<Real> alpha_eff_d;
182 Gpu::AsyncVector<int> eddy_diff_idx_d,eddy_diff_idy_d,eddy_diff_idz_d;
183 alpha_eff_d.resize(alpha_eff.size());
184 eddy_diff_idx_d.resize(eddy_diff_idx.size());
185 eddy_diff_idy_d.resize(eddy_diff_idy.size());
186 eddy_diff_idz_d.resize(eddy_diff_idz.size());
188 Gpu::copy(Gpu::hostToDevice, alpha_eff.begin() , alpha_eff.end() , alpha_eff_d.begin());
189 Gpu::copy(Gpu::hostToDevice, eddy_diff_idx.begin(), eddy_diff_idx.end(), eddy_diff_idx_d.begin());
190 Gpu::copy(Gpu::hostToDevice, eddy_diff_idy.begin(), eddy_diff_idy.end(), eddy_diff_idy_d.begin());
191 Gpu::copy(Gpu::hostToDevice, eddy_diff_idz.begin(), eddy_diff_idz.end(), eddy_diff_idz_d.begin());
194 Real* d_alpha_eff = alpha_eff_d.data();
195 int* d_eddy_diff_idx = eddy_diff_idx_d.data();
196 int* d_eddy_diff_idy = eddy_diff_idy_d.data();
197 int* d_eddy_diff_idz = eddy_diff_idz_d.data();
200 if (l_consA && l_turb) {
201 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
203 const int qty_index = start_comp + n;
204 const int prim_index = qty_index - 1;
206 const int prim_scal_index = prim_index;
218 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
219 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
220 rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_scal_index])
221 + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_scal_index]) );
223 if (ext_dir_on_xlo) {
224 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
225 + 3. * cell_prim(i , j, k, prim_index)
226 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
227 }
else if (ext_dir_on_xhi) {
228 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
229 - 3. * cell_prim(i-1, j, k, prim_index)
230 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
232 xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) -
233 cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
236 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
238 const int qty_index = start_comp + n;
239 const int prim_index = qty_index - 1;
241 const int prim_scal_index = prim_index;
253 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
254 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
255 rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_scal_index])
256 + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_scal_index]) );
258 if (ext_dir_on_ylo) {
259 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
260 + 3. * cell_prim(i, j , k, prim_index)
261 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
262 }
else if (ext_dir_on_yhi) {
263 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
264 - 3. * cell_prim(i, j-1, k, prim_index)
265 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
267 yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
270 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
272 const int qty_index = start_comp + n;
273 const int prim_index = qty_index - 1;
275 const int prim_scal_index = prim_index;
277 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
278 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
279 rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_scal_index])
280 + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_scal_index]) );
291 bool most_on_zlo = ( use_most && exp_most &&
295 if (ext_dir_on_zlo) {
296 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
297 + 3. * cell_prim(i, j, k , prim_index)
298 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
299 }
else if (ext_dir_on_zhi) {
300 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
301 - 3. * cell_prim(i, j, k-1, prim_index)
302 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
304 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
305 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
306 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
308 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
313 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
317 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
320 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
325 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
327 const int qty_index = start_comp + n;
328 const int prim_index = qty_index - 1;
340 Real rhoAlpha = d_alpha_eff[prim_index];
341 rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_index])
342 + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_index]) );
344 if (ext_dir_on_xlo) {
345 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
346 + 3. * cell_prim(i , j, k, prim_index)
347 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
348 }
else if (ext_dir_on_xhi) {
349 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
350 - 3. * cell_prim(i-1, j, k, prim_index)
351 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
353 xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
356 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
358 const int qty_index = start_comp + n;
359 const int prim_index = qty_index - 1;
371 Real rhoAlpha = d_alpha_eff[prim_index];
372 rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_index])
373 + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_index]) );
375 if (ext_dir_on_ylo) {
376 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
377 + 3. * cell_prim(i, j , k, prim_index)
378 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
379 }
else if (ext_dir_on_yhi) {
380 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
381 - 3. * cell_prim(i, j-1, k, prim_index)
382 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
384 yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
387 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
389 const int qty_index = start_comp + n;
390 const int prim_index = qty_index - 1;
392 Real rhoAlpha = d_alpha_eff[prim_index];
393 rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_index])
394 + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_index]) );
405 bool most_on_zlo = ( use_most && exp_most &&
409 if (ext_dir_on_zlo) {
410 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
411 + 3. * cell_prim(i, j, k , prim_index)
412 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
413 }
else if (ext_dir_on_zhi) {
414 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
415 - 3. * cell_prim(i, j, k-1, prim_index)
416 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
418 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
419 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
420 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
422 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
427 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
431 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
434 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
439 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
441 const int qty_index = start_comp + n;
442 const int prim_index = qty_index - 1;
454 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
455 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
457 if (ext_dir_on_xlo) {
458 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
459 + 3. * cell_prim(i , j, k, prim_index)
460 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
461 }
else if (ext_dir_on_xhi) {
462 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
463 - 3. * cell_prim(i-1, j, k, prim_index)
464 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
466 xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
469 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
471 const int qty_index = start_comp + n;
472 const int prim_index = qty_index - 1;
484 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
485 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
487 if (ext_dir_on_ylo) {
488 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
489 + 3. * cell_prim(i, j , k, prim_index)
490 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
491 }
else if (ext_dir_on_yhi) {
492 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
493 - 3. * cell_prim(i, j-1, k, prim_index)
494 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
496 yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
499 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
501 const int qty_index = start_comp + n;
502 const int prim_index = qty_index - 1;
504 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
505 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
516 bool most_on_zlo = ( use_most && exp_most &&
520 if (ext_dir_on_zlo) {
521 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
522 + 3. * cell_prim(i, j, k , prim_index)
523 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
524 }
else if (ext_dir_on_zhi) {
525 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
526 - 3. * cell_prim(i, j, k-1, prim_index)
527 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
529 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
530 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
531 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
533 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
538 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
542 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
545 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
551 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
553 const int qty_index = start_comp + n;
554 const int prim_index = qty_index - 1;
566 Real rhoAlpha = d_alpha_eff[prim_index];
568 if (ext_dir_on_xlo) {
569 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
570 + 3. * cell_prim(i , j, k, prim_index)
571 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
572 }
else if (ext_dir_on_xhi) {
573 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
574 - 3. * cell_prim(i-1, j, k, prim_index)
575 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
577 xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
580 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
582 const int qty_index = start_comp + n;
583 const int prim_index = qty_index - 1;
595 Real rhoAlpha = d_alpha_eff[prim_index];
597 if (ext_dir_on_ylo) {
598 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
599 + 3. * cell_prim(i, j , k, prim_index)
600 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
601 }
else if (ext_dir_on_yhi) {
602 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
603 - 3. * cell_prim(i, j-1, k, prim_index)
604 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
606 yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
609 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
611 const int qty_index = start_comp + n;
612 const int prim_index = qty_index - 1;
614 Real rhoAlpha = d_alpha_eff[prim_index];
625 bool most_on_zlo = ( use_most && exp_most &&
629 if (ext_dir_on_zlo) {
630 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
631 + 3. * cell_prim(i, j, k , prim_index)
632 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
633 }
else if (ext_dir_on_zhi) {
634 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
635 - 3. * cell_prim(i, j, k-1, prim_index)
636 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
638 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
639 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
640 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
642 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
647 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
651 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
654 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
660 for (
int n(0); n < num_comp; ++n)
662 int qty_index = start_comp + n;
663 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
666 cell_rhs(i,j,k,qty_index) -= (xflux(i+1,j ,k ,qty_index) - xflux(i, j, k, qty_index)) * dx_inv * mf_m(i,j,0)
667 +(yflux(i ,j+1,k ,qty_index) - yflux(i, j, k, qty_index)) * dy_inv * mf_m(i,j,0)
668 +(zflux(i ,j ,k+1,qty_index) - zflux(i, j, k, qty_index)) * dz_inv;
687 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
700 cell_rhs(i,j,k,qty_index) += l_abs_g * l_inv_theta0 * hfx_z(i,j,k);
706 cell_rhs(i,j,k,qty_index) += 2.0*mu_turb(i,j,k,
EddyDiff::Mom_v) * SmnSmn_a(i,j,k);
709 cell_rhs(i,j,k,qty_index) -= diss(i,j,k);
714 if (l_use_mynn && start_comp <= RhoKE_comp && end_comp >=
RhoKE_comp) {
718 const int rhoqv_comp = solverChoice.
RhoQv_comp;
719 const int rhoqc_comp = solverChoice.
RhoQc_comp;
720 const int rhoqr_comp = solverChoice.
RhoQr_comp;
722 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
732 cell_rhs(i, j, k, qty_index) +=
ComputeQKESourceTerms(i,j,k,u,v,cell_data,cell_prim,
733 mu_turb,cellSizeInv,domain,
734 pbl_mynn_B1_l,tm_arr(i,j,0),
735 rhoqv_comp, rhoqc_comp, rhoqr_comp,
736 c_ext_dir_on_zlo, c_ext_dir_on_zhi,
737 u_ext_dir_on_zlo, u_ext_dir_on_zhi,
738 v_ext_dir_on_zlo, v_ext_dir_on_zhi,
void DiffusionSrcForState_N(const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &exp_most, 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 GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_m, const Array4< const Real > &mf_u, const Array4< const Real > &mf_v, Array4< Real > &hfx_z, 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_most)
Definition: ERF_DiffusionSrcForState_N.cpp:40
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:53
#define NPRIMVAR_max
Definition: ERF_IndexDefines.H:33
#define PrimQ2_comp
Definition: ERF_IndexDefines.H:54
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimQ4_comp
Definition: ERF_IndexDefines.H:56
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define PrimQ6_comp
Definition: ERF_IndexDefines.H:58
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define PrimQ5_comp
Definition: ERF_IndexDefines.H:57
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define PrimScalar_comp
Definition: ERF_IndexDefines.H:52
#define PrimQ3_comp
Definition: ERF_IndexDefines.H:55
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeQKESourceTerms(int i, int j, int k, const amrex::Array4< const amrex::Real > &uvel, const amrex::Array4< const amrex::Real > &vvel, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &K_turb, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Box &domain, amrex::Real pbl_mynn_B1_l, const amrex::Real theta_mean, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp, bool c_ext_dir_on_zlo, bool c_ext_dir_on_zhi, bool u_ext_dir_on_zlo, bool u_ext_dir_on_zhi, bool v_ext_dir_on_zlo, bool v_ext_dir_on_zhi, const bool use_most=false, const amrex::Real met_h_zeta=1.0)
Definition: ERF_PBLModels.H:160
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:80
@ yvel_bc
Definition: ERF_IndexDefines.H:89
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ xvel_bc
Definition: ERF_IndexDefines.H:88
@ foextrap
Definition: ERF_IndexDefines.H:179
@ ext_dir
Definition: ERF_IndexDefines.H:180
@ ext_dir_prim
Definition: ERF_IndexDefines.H:182
@ Theta_v
Definition: ERF_IndexDefines.H:157
@ Scalar_v
Definition: ERF_IndexDefines.H:159
@ Q_v
Definition: ERF_IndexDefines.H:160
@ Q_h
Definition: ERF_IndexDefines.H:155
@ Mom_v
Definition: ERF_IndexDefines.H:156
@ Scalar_h
Definition: ERF_IndexDefines.H:154
@ KE_v
Definition: ERF_IndexDefines.H:158
@ Theta_h
Definition: ERF_IndexDefines.H:152
@ KE_h
Definition: ERF_IndexDefines.H:153
Definition: ERF_DiffStruct.H:19
amrex::Real rhoAlpha_C
Definition: ERF_DiffStruct.H:95
amrex::Real rhoAlpha_T
Definition: ERF_DiffStruct.H:94
amrex::Real alpha_T
Definition: ERF_DiffStruct.H:87
amrex::Real alpha_C
Definition: ERF_DiffStruct.H:88
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
amrex::Real B1
Definition: ERF_MYNNStruct.H:43
int RhoQr_comp
Definition: ERF_DataStruct.H:700
DiffChoice diffChoice
Definition: ERF_DataStruct.H:592
int RhoQc_comp
Definition: ERF_DataStruct.H:694
int RhoQv_comp
Definition: ERF_DataStruct.H:693
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:594
Definition: ERF_TurbStruct.H:31
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:243
PBLType pbl_type
Definition: ERF_TurbStruct.H:241
RANSType rans_type
Definition: ERF_TurbStruct.H:238
LESType les_type
Definition: ERF_TurbStruct.H:205
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:235