Function for computing the slow RHS for the evolution equations for the density, potential temperature and momentum.
60 BL_PROFILE_REGION(
"erf_make_sources()");
74 const bool l_use_KE = tc.
use_tke;
77 const Box& domain = geom.Domain();
79 const GpuArray<Real, AMREX_SPACEDIM>
dxInv = geom.InvCellSizeArray();
80 const GpuArray<Real, AMREX_SPACEDIM>
dx = geom.CellSizeArray();
94 bool has_moisture = (solverChoice.
moisture_type != MoistureType::None);
99 Table1D<Real> dptr_r_plane, dptr_t_plane, dptr_qv_plane, dptr_qc_plane;
100 TableData<Real, 1> r_plane_tab, t_plane_tab, qv_plane_tab, qc_plane_tab;
103 if (compute_averages)
114 IntVect ng_c(S_data[
IntVars::cons].nGrowVect()); ng_c[2] = 1;
119 int ncomp = (!has_moisture) ? 2 :
RhoQ2_comp+1;
123 cons_ave.compute_averages(
ZDir(), cons_ave.field());
125 int ncell = cons_ave.ncell_line();
127 Gpu::HostVector< Real> r_plane_h(ncell);
128 Gpu::DeviceVector< Real> r_plane_d(ncell);
130 Gpu::HostVector< Real> t_plane_h(ncell);
131 Gpu::DeviceVector< Real> t_plane_d(ncell);
133 cons_ave.line_average(
Rho_comp , r_plane_h);
136 Gpu::copyAsync(Gpu::hostToDevice, r_plane_h.begin(), r_plane_h.end(), r_plane_d.begin());
137 Gpu::copyAsync(Gpu::hostToDevice, t_plane_h.begin(), t_plane_h.end(), t_plane_d.begin());
139 Real* dptr_r = r_plane_d.data();
140 Real* dptr_t = t_plane_d.data();
142 Box tdomain = domain; tdomain.grow(2,ng_c[2]);
143 r_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
144 t_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
148 dptr_r_plane = r_plane_tab.table();
149 dptr_t_plane = t_plane_tab.table();
150 ParallelFor(ncell, [=] AMREX_GPU_DEVICE (
int k) noexcept
152 dptr_r_plane(k-
offset) = dptr_r[k];
153 dptr_t_plane(k-
offset) = dptr_t[k];
158 Gpu::HostVector< Real> qv_plane_h(ncell), qc_plane_h(ncell);
159 Gpu::DeviceVector<Real> qv_plane_d(ncell), qc_plane_d(ncell);
162 cons_ave.line_average(
RhoQ1_comp, qv_plane_h);
163 Gpu::copyAsync(Gpu::hostToDevice, qv_plane_h.begin(), qv_plane_h.end(), qv_plane_d.begin());
166 cons_ave.line_average(
RhoQ2_comp, qc_plane_h);
167 Gpu::copyAsync(Gpu::hostToDevice, qc_plane_h.begin(), qc_plane_h.end(), qc_plane_d.begin());
169 Real* dptr_qv = qv_plane_d.data();
170 Real* dptr_qc = qc_plane_d.data();
172 qv_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
173 qc_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
175 dptr_qv_plane = qv_plane_tab.table();
176 dptr_qc_plane = qc_plane_tab.table();
177 ParallelFor(ncell, [=] AMREX_GPU_DEVICE (
int k) noexcept
179 dptr_qv_plane(k-
offset) = dptr_qv[k];
180 dptr_qc_plane(k-
offset) = dptr_qc[k];
189 int klo = domain.smallEnd(0);
190 int khi = domain.bigEnd(2);
191 int nk =
khi - klo + 2;
192 Gpu::DeviceVector<Real> radiation_flux(nk,
zero);
193 Gpu::DeviceVector<Real> q_integral(nk,
zero);
194 Real* rad_flux = radiation_flux.data();
195 Real* q_int = q_integral.data();
217 #pragma omp parallel if (Gpu::notInLaunchRegion())
222 Box bx = mfi.tilebox();
224 const Array4<const Real>& cell_data = S_data[
IntVars::cons].array(mfi);
225 const Array4<const Real>& cell_prim = S_prim.array(mfi);
226 const Array4<Real> & cell_src = source.array(mfi);
228 const Array4<const Real>&
r0 = r_hse.const_array(mfi);
229 const Array4<const Real>& th0 = th_hse.const_array(mfi);
230 const Array4<const Real>& qv0 = qv_hse.const_array(mfi);
232 const Array4<const Real>& z_cc_arr = z_phys_cc->const_array(mfi);
234 const Array4<const Real>& t_blank_arr = (terrain_blank) ? terrain_blank->const_array(mfi) :
235 Array4<const Real>{};
241 if (solverChoice.
rad_type != RadiationType::None && is_slow_step) {
242 auto const& qheating_arr = qheating_rates->const_array(mfi);
243 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
246 cell_src(i,j,k,
RhoTheta_comp) += cell_data(i,j,k,
Rho_comp) * ( qheating_arr(i,j,k,0) + qheating_arr(i,j,k,1) );
257 if ((is_slow_step && !use_Rayleigh_fast) || (!is_slow_step && use_Rayleigh_fast)) {
262 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
265 Real sinesq = d_sinesq_at_lev[k];
266 cell_src(i, j, k, n) -= dampcoef*sinesq * (
theta -
thetabar[k]) * cell_data(i,j,k,
nr);
276 auto const& rhotheta_src_arr = rhotheta_src->const_array(mfi);
281 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
283 cell_src(i, j, k, n) += cell_data(i,j,k,
nr) * rhotheta_src_arr(i, j, k);
286 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
288 cell_src(i, j, k, n) += rhotheta_src_arr(i, j, k);
294 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
296 cell_src(i, j, k, n) += cell_data(i,j,k,
nr) * rhotheta_src_arr(0, 0, k);
299 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
301 cell_src(i, j, k, n) += rhotheta_src_arr(0, 0, k);
312 auto const& rhoqt_src_arr = rhoqt_src->const_array(mfi);
317 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
319 cell_src(i, j, k, n) += cell_data(i,j,k,
nr) * rhoqt_src_arr(i, j, k);
322 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
324 cell_src(i, j, k, n) += rhoqt_src_arr(i, j, k);
330 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
332 cell_src(i, j, k, n) += cell_data(i,j,k,
nr) * rhoqt_src_arr(0, 0, k);
335 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
337 cell_src(i, j, k, n) += rhoqt_src_arr(0, 0, k);
350 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
352 Real dzInv = (z_cc_arr) ?
one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) :
myhalf*
dxInv[2];
353 Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
354 Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
355 Real wbar_cc =
myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
356 cell_src(i, j, k, n) -= cell_data(i,j,k,
nr) * wbar_cc * (T_hi - T_lo) * dzInv;
359 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
361 Real dzInv = (z_cc_arr) ?
one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) :
myhalf*
dxInv[2];
362 Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
363 Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
364 Real wbar_cc =
myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
365 cell_src(i, j, k, n) -= wbar_cc * (T_hi - T_lo) * dzInv;
377 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
379 Real dzInv = (z_cc_arr) ?
one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) :
myhalf*
dxInv[2];
380 Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
381 Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
382 Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
383 Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
384 Real wbar_cc =
myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
385 cell_src(i, j, k, nv ) -= cell_data(i,j,k,
nr) * wbar_cc * (Qv_hi - Qv_lo) * dzInv;
386 cell_src(i, j, k, nv+1) -= cell_data(i,j,k,
nr) * wbar_cc * (Qc_hi - Qc_lo) * dzInv;
389 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
391 Real dzInv = (z_cc_arr) ?
one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) :
myhalf*
dxInv[2];
392 Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
393 Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
394 Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
395 Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
396 Real wbar_cc =
myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
397 cell_src(i, j, k, nv ) -= wbar_cc * (Qv_hi - Qv_lo) * dzInv;
398 cell_src(i, j, k, nv+1) -= wbar_cc * (Qc_hi - Qc_lo) * dzInv;
406 if (l_use_ndiff && is_slow_step)
408 const Array4<const Real>& mf_mx = mapfac[
MapFacType::m_x]->const_array(mfi);
409 const Array4<const Real>& mf_my = mapfac[
MapFacType::m_y]->const_array(mfi);
413 cell_data, cell_data, cell_src, mf_mx, mf_my);
417 cell_prim, cell_data, cell_src, mf_mx, mf_my);
420 if (l_use_KE && l_diff_KE) {
422 cell_prim, cell_data, cell_src, mf_mx, mf_my);
426 cell_prim, cell_data, cell_src, mf_mx, mf_my);
438 const Array4<const Real>& surface_state_arr = (*surface_state_at_lev).array(mfi);
448 const amrex::Array4<const amrex::Real>& pert_cell = turbPert.
pb_cell[level].const_array(mfi);
466 for (
int nt = 1; nt < n_sounding_times; nt++) {
469 if (itime_n == n_sounding_times-1) {
472 itime_np1 = itime_n+1;
475 coeff_n =
one - coeff_np1;
484 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
486 Real nudge = (coeff_n*theta_inp_sound_n[k] + coeff_np1*theta_inp_sound_np1[k]) - (dptr_t_plane(k)/dptr_r_plane(k));
488 cell_src(i, j, k, n) += cell_data(i, j, k,
nr) * nudge;
495 if (solverChoice.
terrain_type == TerrainType::ImmersedForcing &&
496 ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
498 const Array4<const Real>& u =
xvel.array(mfi);
499 const Array4<const Real>& v =
yvel.array(mfi);
502 const Real* dx_arr = geom.CellSize();
503 const Real dx_x = dx_arr[0];
504 const Real dx_y = dx_arr[1];
505 const Real dx_z = dx_arr[2];
508 const Real drag_coefficient = alpha_h / std::pow(dx_x*dx_y*dx_z,
one/
three);
525 ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
527 const Real t_blank = t_blank_arr(i, j, k);
528 const Real t_blank_above = t_blank_arr(i, j, k+1);
529 const Real ux_cc_2r =
myhalf * (u(i ,j ,k+1) + u(i+1,j ,k+1));
530 const Real uy_cc_2r =
myhalf * (v(i ,j ,k+1) + v(i ,j+1,k+1));
531 const Real h_windspeed2r = std::sqrt(ux_cc_2r * ux_cc_2r + uy_cc_2r * uy_cc_2r);
537 if (init_surf_temp >
zero) {
538 if (t_blank > 0 && (t_blank_above ==
zero)) {
539 const Real surf_temp = init_surf_temp + surf_heating_rate*time;
541 cell_src(i, j, k-1,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
546 if (tflux !=
Real(1e-8)){
547 if (t_blank > 0 && (t_blank_above ==
zero)) {
551 Real ustar = h_windspeed2r * kappa / (std::log((
Real(1.5)) * dx_z /
z0) - psi_m);
552 const Real Olen = -ustar * ustar * ustar *
theta / (kappa * ggg * tflux + tiny);
554 const Real zeta_neighbor = (
Real(1.5)) * dx_z / Olen;
559 psi_h_neighbor = sfuns.
calc_psi_h(zeta_neighbor);
560 ustar = h_windspeed2r * kappa / (std::log((
Real(1.5)) * dx_z /
z0) - psi_m);
563 if (!(ustar >
zero && !std::isnan(ustar))) { ustar =
zero; }
564 if (!(ustar <
two && !std::isnan(ustar))) { ustar =
two; }
565 if (psi_h_neighbor > std::log(
Real(1.5) * dx_z /
z0)) { psi_h_neighbor = std::log(
Real(1.5) * dx_z /
z0); }
566 if (psi_h > std::log(
myhalf * dx_z /
z0)) { psi_h = std::log(
myhalf * dx_z /
z0); }
569 const Real thetastar =
theta * ustar * ustar / (kappa * ggg * Olen);
570 const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((
Real(1.5)) * dx_z /
z0) - psi_h_neighbor);
571 const Real tTarget = surf_temp + thetastar / kappa * (std::log((
myhalf) * dx_z /
z0) - psi_h);
574 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
579 if (Olen_in !=
Real(1e-8)){
580 if (t_blank > 0 && (t_blank_above ==
zero)) {
581 const Real Olen = Olen_in;
583 const Real zeta_neighbor = (
Real(1.5)) * dx_z / Olen;
589 const Real ustar = h_windspeed2r * kappa / (std::log((
Real(1.5)) * dx_z /
z0) - psi_m);
592 const Real thetastar =
theta * ustar * ustar / (kappa * ggg * Olen);
593 const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((
Real(1.5)) * dx_z /
z0) - psi_h_neighbor);
594 const Real tTarget = surf_temp + thetastar / kappa * (std::log((
myhalf) * dx_z /
z0) - psi_h);
597 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
607 if ((solverChoice.
buildings_type == BuildingsType::ImmersedForcing ) &&
608 ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
611 const Real* dx_arr = geom.CellSize();
612 const Real dx_x = dx_arr[0];
613 const Real dx_y = dx_arr[1];
617 const Real min_t_blank =
Real(0.005);
624 ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) noexcept
626 Real t_blank = t_blank_arr(i, j, k);
627 Real t_blank_below = t_blank_arr(i, j, k-1);
628 Real t_blank_above = t_blank_arr(i, j, k+1);
629 Real t_blank_north = t_blank_arr(i , j+1, k);
630 Real t_blank_south = t_blank_arr(i , j-1, k);
631 Real t_blank_east = t_blank_arr(i+1, j , k);
632 Real t_blank_west = t_blank_arr(i-1, j , k);
633 if (t_blank < min_t_blank) { t_blank =
zero; }
634 if (t_blank_below < min_t_blank) { t_blank_below =
zero; }
635 if (t_blank_north < min_t_blank) { t_blank_north =
zero; }
636 if (t_blank_south < min_t_blank) { t_blank_south =
zero; }
637 if (t_blank_east < min_t_blank) { t_blank_east =
zero; }
638 if (t_blank_west < min_t_blank) { t_blank_west =
zero; }
640 Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) :
dx[2];
641 Real drag_coefficient = alpha_h / std::pow(dx_x*dx_y*dx_z,
one/
three);
644 if (init_surf_temp >
zero) {
645 const Real surf_temp = init_surf_temp + surf_heating_rate*time;
646 if (t_blank > 0 && (t_blank_above ==
zero) && (t_blank_below ==
one)) {
648 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
650 }
else if (((t_blank > 0 && t_blank < t_blank_west && t_blank_east ==
zero) ||
651 (t_blank > 0 && t_blank < t_blank_east && t_blank_west ==
zero) ||
652 (t_blank > 0 && t_blank < t_blank_north && t_blank_south ==
zero) ||
653 (t_blank > 0 && t_blank < t_blank_south && t_blank_north ==
zero))) {
658 if ((t_blank < t_blank_north) && (t_blank_north ==
one)) {
660 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
664 if ((t_blank < t_blank_south) && (t_blank_south ==
one)) {
666 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
670 if ((t_blank < t_blank_east) && (t_blank_east ==
one)) {
672 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
676 if ((t_blank < t_blank_west) && (t_blank_west ==
one)) {
678 cell_src(i, j, k,
RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
698 Box xybx = makeSlab(bx,2,klo);
699 ParallelFor(xybx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int ) noexcept
705 for (
int k(klo+1); k<=
khi+1; ++k) {
708 Real dz = (z_cc_arr) ?
myhalf * (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-2)) :
dx[2];
709 q_int[lk] = q_int[lk-1] + krad * cell_data(i,j,k-1,
Rho_comp) * cell_data(i,j,k-1,
RhoQ2_comp) *
dz;
712 if ( (qt_lo > qt_i) && (qt_hi < qt_i) ) {
713 zi =
myhalf * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
720 for (
int k(klo); k<=
khi+1; ++k) {
722 Real z =
myhalf * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
723 rad_flux[lk] = F1*std::exp(-q_int[lk]) + F0*std::exp(-(q_int_inf - q_int[lk]));
730 for (
int k(klo); k<=
khi; ++k) {
733 Real dzInv = (z_cc_arr) ?
one/ (
myhalf * (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1))) :
dxInv[2];
736 Real dTdt = (rad_flux[lk+1] - rad_flux[lk]) * dzInv / (-cell_data(i,j,k,
Rho_comp)*
Cp_d);
void ApplySpongeZoneBCsForCC(const SpongeChoice &spongeChoice, const Geometry geom, const Box &bx, const Array4< Real > &cell_rhs, const Array4< const Real > &cell_data, const Array4< const Real > &r0, const Array4< const Real > &th0, const Array4< const Real > &qv0, const Array4< const Real > &z_phys_cc, int n_qstate)
Definition: ERF_ApplySpongeZoneBCs.cpp:7
void ApplySurfaceTreatment_BulkCoeff_CC(const Box &bx, const Array4< Real > &cell_rhs, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_cc, const Array4< const Real > &surface_state_arr)
Definition: ERF_ApplySurfaceTreatment_BulkCoeff.cpp:60
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:39
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:31
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:40
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
@ thetabar
Definition: ERF_DataStruct.H:98
@ m_y
Definition: ERF_DataStruct.H:25
@ m_x
Definition: ERF_DataStruct.H:24
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:156
#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 NDRY
Definition: ERF_IndexDefines.H:13
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:55
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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);})
void NumericalDiffusion_Scal(const Box &bx, const int start_comp, const int num_comp, const Real dt, const Real num_diff_coeff, const Array4< const Real > &prim_data, const Array4< const Real > &cell_data, const Array4< Real > &rhs, const Array4< const Real > &mfx_arr, const Array4< const Real > &mfy_arr)
Definition: ERF_NumericalDiffusion.cpp:18
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
Definition: ERF_PlaneAverage.H:14
@ qv0_comp
Definition: ERF_IndexDefines.H:77
@ th0_comp
Definition: ERF_IndexDefines.H:76
@ r0_comp
Definition: ERF_IndexDefines.H:73
@ cons
Definition: ERF_IndexDefines.H:192
@ theta
Definition: ERF_MM5.H:20
@ qv
Definition: ERF_Kessler.H:29
@ nr
Definition: ERF_Morrison.H:45
@ xvel
Definition: ERF_IndexDefines.H:175
@ cons
Definition: ERF_IndexDefines.H:174
@ yvel
Definition: ERF_IndexDefines.H:176
@ dz
Definition: ERF_AdvanceWSM6.cpp:104
@ zi
Definition: ERF_AdvanceWSM6.cpp:133
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
real(c_double), private rhoi
Definition: ERF_module_mp_morr_two_moment.F90:188
real(kind=kind_phys), parameter, private r0
Definition: ERF_module_mp_wsm6.F90:21
bool rayleigh_damp_T
Definition: ERF_DampingStruct.H:87
amrex::Real rayleigh_dampcoef
Definition: ERF_DampingStruct.H:88
RayleighDampingType rayleigh_damping_type
Definition: ERF_DampingStruct.H:101
amrex::Real if_init_surf_temp
Definition: ERF_DataStruct.H:1208
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
amrex::Real if_surf_temp_flux
Definition: ERF_DataStruct.H:1207
DampingChoice dampingChoice
Definition: ERF_DataStruct.H:1144
bool spatial_moisture_forcing
Definition: ERF_DataStruct.H:1243
RadiationType rad_type
Definition: ERF_DataStruct.H:1303
amrex::Real if_z0
Definition: ERF_DataStruct.H:1206
bool do_theta_advection
Definition: ERF_DataStruct.H:1238
bool custom_rhotheta_forcing
Definition: ERF_DataStruct.H:1235
bool custom_w_subsidence
Definition: ERF_DataStruct.H:1237
bool nudging_from_input_sounding
Definition: ERF_DataStruct.H:1247
amrex::Real if_Cd_scalar
Definition: ERF_DataStruct.H:1203
bool immersed_forcing_substep
Definition: ERF_DataStruct.H:1199
amrex::Real if_Olen_in
Definition: ERF_DataStruct.H:1210
bool use_num_diff
Definition: ERF_DataStruct.H:1295
bool four_stream_radiation
Definition: ERF_DataStruct.H:1196
bool custom_moisture_forcing
Definition: ERF_DataStruct.H:1236
amrex::Real num_diff_coeff
Definition: ERF_DataStruct.H:1296
amrex::Real if_surf_heating_rate
Definition: ERF_DataStruct.H:1209
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1146
MoistureType moisture_type
Definition: ERF_DataStruct.H:1299
bool custom_forcing_prim_vars
Definition: ERF_DataStruct.H:1241
static BuildingsType buildings_type
Definition: ERF_DataStruct.H:1128
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1125
SpongeChoice spongeChoice
Definition: ERF_DataStruct.H:1145
static InitType init_type
Definition: ERF_DataStruct.H:1119
bool use_source_perturbation(int lev) const
Definition: ERF_DataStruct.H:1274
bool hindcast_surface_bcs
Definition: ERF_DataStruct.H:1345
bool spatial_rhotheta_forcing
Definition: ERF_DataStruct.H:1242
int ave_plane
Definition: ERF_DataStruct.H:1313
static SpongeType sponge_type
Definition: ERF_SpongeStruct.H:81
Definition: ERF_TurbStruct.H:82
bool use_tke
Definition: ERF_TurbStruct.H:471
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:506
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:667
void apply_tpi(const int &lev, const amrex::Box &vbx, const int &comp, const amrex::IndexType &m_ixtype, const amrex::Array4< amrex::Real > &src_arr, const amrex::Array4< amrex::Real const > &pert_cell)
Definition: ERF_TurbPertStruct.H:351
Definition: ERF_MOSTStress.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_m(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:90
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_h(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:104