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_ddorf = (turbChoice.
les_type == LESType::Deardorff);
86 bool l_use_mynn = (turbChoice.
pbl_type == PBLType::MYNN25);
89 bool l_turb = ( (turbChoice.
les_type == LESType::Smagorinsky) ||
90 (turbChoice.
les_type == LESType::Deardorff ) ||
91 (turbChoice.
pbl_type == PBLType::MYNN25 ) ||
92 (turbChoice.
pbl_type == PBLType::YSU ) );
94 const Box xbx = surroundingNodes(bx,0);
95 const Box ybx = surroundingNodes(bx,1);
96 const Box zbx = surroundingNodes(bx,2);
98 const int end_comp = start_comp + num_comp - 1;
179 Gpu::AsyncVector<Real> alpha_eff_d;
180 Gpu::AsyncVector<int> eddy_diff_idx_d,eddy_diff_idy_d,eddy_diff_idz_d;
181 alpha_eff_d.resize(alpha_eff.size());
182 eddy_diff_idx_d.resize(eddy_diff_idx.size());
183 eddy_diff_idy_d.resize(eddy_diff_idy.size());
184 eddy_diff_idz_d.resize(eddy_diff_idz.size());
186 Gpu::copy(Gpu::hostToDevice, alpha_eff.begin() , alpha_eff.end() , alpha_eff_d.begin());
187 Gpu::copy(Gpu::hostToDevice, eddy_diff_idx.begin(), eddy_diff_idx.end(), eddy_diff_idx_d.begin());
188 Gpu::copy(Gpu::hostToDevice, eddy_diff_idy.begin(), eddy_diff_idy.end(), eddy_diff_idy_d.begin());
189 Gpu::copy(Gpu::hostToDevice, eddy_diff_idz.begin(), eddy_diff_idz.end(), eddy_diff_idz_d.begin());
192 Real* d_alpha_eff = alpha_eff_d.data();
193 int* d_eddy_diff_idx = eddy_diff_idx_d.data();
194 int* d_eddy_diff_idy = eddy_diff_idy_d.data();
195 int* d_eddy_diff_idz = eddy_diff_idz_d.data();
198 if (l_consA && l_turb) {
199 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
201 const int qty_index = start_comp + n;
202 const int prim_index = qty_index - 1;
204 const int prim_scal_index = prim_index;
216 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
217 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
218 rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_scal_index])
219 + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_scal_index]) );
221 if (ext_dir_on_xlo) {
222 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
223 + 3. * cell_prim(i , j, k, prim_index)
224 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
225 }
else if (ext_dir_on_xhi) {
226 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
227 - 3. * cell_prim(i-1, j, k, prim_index)
228 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
230 xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) -
231 cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
234 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
236 const int qty_index = start_comp + n;
237 const int prim_index = qty_index - 1;
239 const int prim_scal_index = prim_index;
251 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
252 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
253 rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_scal_index])
254 + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_scal_index]) );
256 if (ext_dir_on_ylo) {
257 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
258 + 3. * cell_prim(i, j , k, prim_index)
259 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
260 }
else if (ext_dir_on_yhi) {
261 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
262 - 3. * cell_prim(i, j-1, k, prim_index)
263 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
265 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);
268 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
270 const int qty_index = start_comp + n;
271 const int prim_index = qty_index - 1;
273 const int prim_scal_index = prim_index;
275 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
276 Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
277 rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_scal_index])
278 + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_scal_index]) );
289 bool most_on_zlo = ( use_most && exp_most &&
293 if (ext_dir_on_zlo) {
294 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
295 + 3. * cell_prim(i, j, k , prim_index)
296 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
297 }
else if (ext_dir_on_zhi) {
298 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
299 - 3. * cell_prim(i, j, k-1, prim_index)
300 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
302 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
303 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
304 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
306 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
311 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
315 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
318 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
323 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
325 const int qty_index = start_comp + n;
326 const int prim_index = qty_index - 1;
338 Real rhoAlpha = d_alpha_eff[prim_index];
339 rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_index])
340 + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_index]) );
342 if (ext_dir_on_xlo) {
343 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
344 + 3. * cell_prim(i , j, k, prim_index)
345 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
346 }
else if (ext_dir_on_xhi) {
347 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
348 - 3. * cell_prim(i-1, j, k, prim_index)
349 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
351 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);
354 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
356 const int qty_index = start_comp + n;
357 const int prim_index = qty_index - 1;
369 Real rhoAlpha = d_alpha_eff[prim_index];
370 rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_index])
371 + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_index]) );
373 if (ext_dir_on_ylo) {
374 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
375 + 3. * cell_prim(i, j , k, prim_index)
376 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
377 }
else if (ext_dir_on_yhi) {
378 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
379 - 3. * cell_prim(i, j-1, k, prim_index)
380 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
382 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);
385 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
387 const int qty_index = start_comp + n;
388 const int prim_index = qty_index - 1;
390 Real rhoAlpha = d_alpha_eff[prim_index];
391 rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_index])
392 + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_index]) );
403 bool most_on_zlo = ( use_most && exp_most &&
407 if (ext_dir_on_zlo) {
408 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
409 + 3. * cell_prim(i, j, k , prim_index)
410 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
411 }
else if (ext_dir_on_zhi) {
412 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
413 - 3. * cell_prim(i, j, k-1, prim_index)
414 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
416 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
417 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
418 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
420 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
425 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
429 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
432 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
437 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
439 const int qty_index = start_comp + n;
440 const int prim_index = qty_index - 1;
452 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i-1, j, k,
Rho_comp) );
453 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
455 if (ext_dir_on_xlo) {
456 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
457 + 3. * cell_prim(i , j, k, prim_index)
458 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
459 }
else if (ext_dir_on_xhi) {
460 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
461 - 3. * cell_prim(i-1, j, k, prim_index)
462 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
464 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);
467 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
469 const int qty_index = start_comp + n;
470 const int prim_index = qty_index - 1;
482 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j-1, k,
Rho_comp) );
483 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
485 if (ext_dir_on_ylo) {
486 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
487 + 3. * cell_prim(i, j , k, prim_index)
488 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
489 }
else if (ext_dir_on_yhi) {
490 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
491 - 3. * cell_prim(i, j-1, k, prim_index)
492 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
494 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);
497 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
499 const int qty_index = start_comp + n;
500 const int prim_index = qty_index - 1;
502 Real rhoFace = 0.5 * ( cell_data(i, j, k,
Rho_comp) + cell_data(i, j, k-1,
Rho_comp) );
503 Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
514 bool most_on_zlo = ( use_most && exp_most &&
518 if (ext_dir_on_zlo) {
519 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
520 + 3. * cell_prim(i, j, k , prim_index)
521 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
522 }
else if (ext_dir_on_zhi) {
523 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
524 - 3. * cell_prim(i, j, k-1, prim_index)
525 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
527 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
528 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
529 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
531 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
536 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
540 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
543 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
549 ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
551 const int qty_index = start_comp + n;
552 const int prim_index = qty_index - 1;
564 Real rhoAlpha = d_alpha_eff[prim_index];
566 if (ext_dir_on_xlo) {
567 xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
568 + 3. * cell_prim(i , j, k, prim_index)
569 - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
570 }
else if (ext_dir_on_xhi) {
571 xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
572 - 3. * cell_prim(i-1, j, k, prim_index)
573 + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
575 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);
578 ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
580 const int qty_index = start_comp + n;
581 const int prim_index = qty_index - 1;
593 Real rhoAlpha = d_alpha_eff[prim_index];
595 if (ext_dir_on_ylo) {
596 yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
597 + 3. * cell_prim(i, j , k, prim_index)
598 - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
599 }
else if (ext_dir_on_yhi) {
600 yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
601 - 3. * cell_prim(i, j-1, k, prim_index)
602 + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
604 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);
607 ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
609 const int qty_index = start_comp + n;
610 const int prim_index = qty_index - 1;
612 Real rhoAlpha = d_alpha_eff[prim_index];
623 bool most_on_zlo = ( use_most && exp_most &&
627 if (ext_dir_on_zlo) {
628 zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
629 + 3. * cell_prim(i, j, k , prim_index)
630 - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
631 }
else if (ext_dir_on_zhi) {
632 zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
633 - 3. * cell_prim(i, j, k-1, prim_index)
634 + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
636 zflux(i,j,k,qty_index) = hfx_z(i,j,0);
637 }
else if (most_on_zlo && (qty_index ==
RhoQ1_comp)) {
638 zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
640 zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
645 hfx_z(i,j,k) = zflux(i,j,k,qty_index);
649 qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
652 qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
658 for (
int n(0); n < num_comp; ++n)
660 int qty_index = start_comp + n;
661 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
664 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)
665 +(yflux(i ,j+1,k ,qty_index) - yflux(i, j, k, qty_index)) * dy_inv * mf_m(i,j,0)
666 +(zflux(i ,j ,k+1,qty_index) - zflux(i, j, k, qty_index)) * dz_inv;
684 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
697 cell_rhs(i,j,k,qty_index) += l_abs_g * l_inv_theta0 * hfx_z(i,j,k);
703 cell_rhs(i,j,k,qty_index) += 2.0*mu_turb(i,j,k,
EddyDiff::Mom_v) * SmnSmn_a(i,j,k);
706 cell_rhs(i,j,k,qty_index) -= diss(i,j,k);
711 if (l_use_mynn && start_comp <= RhoKE_comp && end_comp >=
RhoKE_comp) {
715 const int rhoqv_comp = solverChoice.
RhoQv_comp;
716 const int rhoqc_comp = solverChoice.
RhoQc_comp;
717 const int rhoqr_comp = solverChoice.
RhoQr_comp;
719 ParallelFor(bx,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
729 cell_rhs(i, j, k, qty_index) +=
ComputeQKESourceTerms(i,j,k,u,v,cell_data,cell_prim,
730 mu_turb,cellSizeInv,domain,
731 pbl_mynn_B1_l,tm_arr(i,j,0),
732 rhoqv_comp, rhoqc_comp, rhoqr_comp,
733 c_ext_dir_on_zlo, c_ext_dir_on_zhi,
734 u_ext_dir_on_zlo, u_ext_dir_on_zhi,
735 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:687
DiffChoice diffChoice
Definition: ERF_DataStruct.H:579
int RhoQc_comp
Definition: ERF_DataStruct.H:681
int RhoQv_comp
Definition: ERF_DataStruct.H:680
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:581
Definition: ERF_TurbStruct.H:31
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:203
PBLType pbl_type
Definition: ERF_TurbStruct.H:201
LESType les_type
Definition: ERF_TurbStruct.H:175
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:195