7 Vector<MultiFab>& S_old,
8 Vector<MultiFab>& S_data,
9 const Real old_step_time,
10 const Real old_stage_time,
11 const Real new_stage_time,
25 BL_PROFILE(
"slow_rhs_fun_pre");
26 if (verbose) Print() << std::setprecision(timeprecision)
27 <<
"Making slow rhs at time " << old_stage_time
28 <<
" for fast variables advancing from " << old_step_time
29 <<
" to " << new_stage_time << std::endl;
31 Real slow_dt = new_stage_time - old_step_time;
36 YAFluxRegister* fr_as_crse =
nullptr;
37 YAFluxRegister* fr_as_fine =
nullptr;
38 if (solverChoice.coupling_type == CouplingType::TwoWay && finest_level > 0) {
39 if (level < finest_level) {
40 fr_as_crse = getAdvFluxReg(level+1);
44 fr_as_fine = getAdvFluxReg(level);
49 MultiFab* forest_drag = (solverChoice.do_forest_drag) ?
50 m_forest_drag[level]->get_drag_field() :
nullptr;
53 MultiFab* terrain_blank = (solverChoice.terrain_type == TerrainType::ImmersedForcing) ?
54 terrain_blanking[level].get() :
nullptr;
58 if (solverChoice.moisture_type != MoistureType::None) {
62 MultiFab p0_to_use, base_to_use;
63 MultiFab *zpn_to_use, *zpc_to_use, *ax_to_use, *ay_to_use, *az_to_use, *dJ_to_use;
66 std::unique_ptr<MultiFab> z_t_pert;
67 if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
74 zpn_to_use = z_phys_nd_src[level].get();
75 zpc_to_use = z_phys_cc_src[level].get();
76 ax_to_use = ax_src[level].get();
77 ay_to_use = ay_src[level].get();
78 az_to_use = az_src[level].get();
79 dJ_to_use = detJ_cc_src[level].get();
84 zpn_to_use = z_phys_nd[level].get();
85 zpc_to_use = z_phys_cc[level].get();
86 ax_to_use = ax[level].get();
87 ay_to_use = ay[level].get();
88 az_to_use = az[level].get();
89 dJ_to_use = detJ_cc[level].get();
96 S_data, S_prim, cc_src, base_state[level], zpc_to_use,
97 qheating_rates[level].get(),
98 terrain_blank, fine_geom, solverChoice,
100 dptr_rhotheta_src, dptr_rhoqt_src,
101 dptr_wbar_sub, d_rayleigh_ptrs_at_lev,
102 input_sounding_data, turbPert,
true);
108 p0_to_use, *zpn_to_use, *zpc_to_use, get_eb(level), gradp[level]);
113 make_buoyancy(S_data, S_prim,
qt, buoyancy, fine_geom, solverChoice, base_to_use,
114 micro->Get_Qstate_Moist_Size(), get_eb(level), solverChoice.anelastic[level]);
120 S_data, *zpn_to_use, *zpc_to_use, stretched_dz_h[level],
121 xvel_new, yvel_new, zvel_new,
122 xmom_src, ymom_src, zmom_src,
123 base_to_use, forest_drag, terrain_blank,
124 cosPhi_m[level].get(), sinPhi_m[level].get(), fine_geom, solverChoice,
126 (solverChoice.have_geo_wind_profile) ? d_u_geos[level].data():
nullptr,
127 (solverChoice.have_geo_wind_profile) ? d_v_geos[level].data():
nullptr,
128 dptr_wbar_sub, d_rayleigh_ptrs_at_lev, d_sponge_ptrs_at_lev,
129 (solverChoice.hindcast_lateral_forcing? &forecast_state_interp[level] :
nullptr),
130 input_sounding_data,
true);
136 xflux_imask[level], yflux_imask[level], zflux_imask[level],
137 thin_xforce[level], thin_yforce[level], thin_zforce[level]);
143 S_prim,
qt, avg_xmom[level], avg_ymom[level], avg_zmom[level],
144 xvel_new, yvel_new, zvel_new,
145 z_t_rk[level], cc_src, xmom_src, ymom_src, zmom_src, buoyancy,
146 (level > 0) ? &zmom_crse_rhs[level] :
nullptr,
147 Tau[level], SmnSmn, eddyDiffs, Hfx1, Hfx2, Hfx3, Q1fx1, Q1fx2, Q1fx3, Q2fx3, Diss,
148 fine_geom, solverChoice, m_SurfaceLayer, domain_bcs_type_d, domain_bcs_type,
149 *zpn_to_use, *zpc_to_use, *ax_to_use, *ay_to_use, *az_to_use, *dJ_to_use,
150 stretched_dz_d[level], gradp[level],
151 mapfac[level], get_eb(level), fr_as_crse, fr_as_fine);
156 if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
163 MultiFab* r0_new = &r_hse_new;
164 MultiFab* p0_new = &p_hse_new;
165 MultiFab* pi0_new = &pi_hse_new;
166 MultiFab* th0_new = &th_hse_new;
171 MultiFab rt0(
p0->boxArray(),
p0->DistributionMap(),1,1);
172 MultiFab rt0_new(
p0->boxArray(),
p0->DistributionMap(),1,1);
173 MultiFab r0_temp(
p0->boxArray(),
p0->DistributionMap(),1,1);
180 Real dt_base = (new_stage_time - old_step_time);
182 const GpuArray<Real, AMREX_SPACEDIM> dxInv = fine_geom.InvCellSizeArray();
184 const Real l_rdOcp = solverChoice.rdOcp;
187 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
189 for ( MFIter mfi(*
p0,TilingIfNotGPU()); mfi.isValid(); ++mfi)
191 const Array4<Real > rt0_arr = rt0.array(mfi);
192 const Array4<Real > rt0_tmp_arr = rt0_new.array(mfi);
194 const Array4<Real const> r0_arr = r0->const_array(mfi);
195 const Array4<Real > r0_new_arr = r0_new->array(mfi);
196 const Array4<Real > r0_tmp_arr = r0_temp.array(mfi);
198 const Array4<Real const> p0_arr =
p0->const_array(mfi);
199 const Array4<Real > p0_new_arr = p0_new->array(mfi);
200 const Array4<Real > pi0_new_arr = pi0_new->array(mfi);
201 const Array4<Real > th0_new_arr = th0_new->array(mfi);
203 const Array4<Real >& z_t_arr = z_t_rk[level]->array(mfi);
205 const Array4<Real const>& dJ_old_arr = detJ_cc[level]->const_array(mfi);
206 const Array4<Real const>& dJ_new_arr = detJ_cc_new[level]->const_array(mfi);
207 const Array4<Real const>& dJ_src_arr = detJ_cc_src[level]->const_array(mfi);
209 Box gbx = mfi.growntilebox({1,1,1});
210 amrex::ParallelFor(gbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
214 r0_tmp_arr(i,j,k) = r0_new_arr(i,j,k);
217 Box gbx2 = mfi.growntilebox({1,1,0});
218 amrex::ParallelFor(gbx2, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
220 Real zflux_r_lo = -z_t_arr(i,j,k ) * 0.5 * (r0_tmp_arr(i,j,k) + r0_tmp_arr(i,j,k-1));
221 Real zflux_r_hi = -z_t_arr(i,j,k+1) * 0.5 * (r0_tmp_arr(i,j,k) + r0_tmp_arr(i,j,k+1));
223 Real zflux_rt_lo = zflux_r_lo * 0.5 * (rt0_tmp_arr(i,j,k)/r0_tmp_arr(i,j,k) + rt0_tmp_arr(i,j,k-1)/r0_tmp_arr(i,j,k-1));
224 Real zflux_rt_hi = zflux_r_hi * 0.5 * (rt0_tmp_arr(i,j,k)/r0_tmp_arr(i,j,k) + rt0_tmp_arr(i,j,k+1)/r0_tmp_arr(i,j,k+1));
226 Real invdetJ = 1.0 / dJ_src_arr(i,j,k);
228 Real src_r = - invdetJ * ( zflux_r_hi - zflux_r_lo ) * dxInv[2];
229 Real src_rt = - invdetJ * ( zflux_rt_hi - zflux_rt_lo ) * dxInv[2];
231 Real rho0_new = dJ_old_arr(i,j,k) * r0_arr(i,j,k) + dt_base * dJ_src_arr(i,j,k) * src_r;
232 Real rt0_tmp_new = dJ_old_arr(i,j,k) * rt0_arr(i,j,k) + dt_base * dJ_src_arr(i,j,k) * src_rt;
234 r0_new_arr(i,j,k) = rho0_new / dJ_new_arr(i,j,k);
235 rt0_tmp_new /= dJ_new_arr(i,j,k);
239 th0_new_arr(i,j,k) = rt0_tmp_new / r0_new_arr(i,j,k);
242 r0_new->FillBoundary(fine_geom.periodicity());
243 p0_new->FillBoundary(fine_geom.periodicity());
244 th0_new->FillBoundary(fine_geom.periodicity());
247 #ifdef ERF_USE_NETCDF
249 if (solverChoice.use_real_bcs && (level == 0)) {
252 new_stage_time, slow_dt, stop_time-start_time,
253 real_width, real_set_width, fine_geom,
254 S_rhs, S_old, S_data,
255 bdy_data_xlo, bdy_data_xhi,
256 bdy_data_ylo, bdy_data_yhi);
void add_thin_body_sources(MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, std::unique_ptr< iMultiFab > &xflux_imask_lev, std::unique_ptr< iMultiFab > &yflux_imask_lev, std::unique_ptr< iMultiFab > &zflux_imask_lev, std::unique_ptr< MultiFab > &thin_xforce_lev, std::unique_ptr< MultiFab > &thin_yforce_lev, std::unique_ptr< MultiFab > &thin_zforce_lev)
Definition: ERF_AddThinBodySources.cpp:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:156
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhoThetagivenP(const amrex::Real p, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:172
void realbdy_compute_interior_ghost_rhs(const Real &bdy_time_interval, const Real &time, const Real &delta_t, const Real &stop_time_elapsed, int width, int set_width, const Geometry &geom, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old_data, Vector< MultiFab > &S_cur_data, Vector< Vector< FArrayBox >> &bdy_data_xlo, Vector< Vector< FArrayBox >> &bdy_data_xhi, Vector< Vector< FArrayBox >> &bdy_data_ylo, Vector< Vector< FArrayBox >> &bdy_data_yhi)
Definition: ERF_InteriorGhostCells.cpp:174
void make_buoyancy(const Vector< MultiFab > &S_data, const MultiFab &S_prim, const MultiFab &qt, MultiFab &buoyancy, const Geometry geom, const SolverChoice &solverChoice, const MultiFab &base_state, const int n_qstate, const eb_ &ebfact, const int anelastic)
Definition: ERF_MakeBuoyancy.cpp:32
void make_gradp_pert(int level, const SolverChoice &solverChoice, const Geometry &geom, Vector< MultiFab > &S_data, MultiFab &p0, MultiFab &z_phys_nd, MultiFab &z_phys_cc, const eb_ &ebfact, Vector< MultiFab > &gradp)
Definition: ERF_MakeGradP.cpp:28
void make_mom_sources(Real time, Real dt, const Vector< MultiFab > &S_data, MultiFab &z_phys_nd, MultiFab &z_phys_cc, Vector< Real > &stretched_dz_h, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &wvel, MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, const MultiFab &base_state, MultiFab *forest_drag, MultiFab *terrain_blank, MultiFab *cosPhi_mf, MultiFab *sinPhi_mf, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &, const Real *dptr_u_geos, const Real *dptr_v_geos, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, const Vector< Real * > d_sponge_ptrs_at_lev, const Vector< MultiFab > *forecast_state_at_lev, InputSoundingData &input_sounding_data, bool is_slow_step)
Definition: ERF_MakeMomSources.cpp:35
void make_sources(int level, int, Real dt, Real time, const Vector< MultiFab > &S_data, const MultiFab &S_prim, MultiFab &source, const MultiFab &base_state, const MultiFab *z_phys_cc, const MultiFab *qheating_rates, MultiFab *terrain_blank, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &mapfac, const Real *dptr_rhotheta_src, const Real *dptr_rhoqt_src, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, InputSoundingData &input_sounding_data, TurbulentPerturbation &turbPert, bool is_slow_step)
Definition: ERF_MakeSources.cpp:31
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void erf_slow_rhs_pre(int level, int finest_level, int nrk, Real dt, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old, Vector< MultiFab > &S_data, const MultiFab &S_prim, const MultiFab &qt, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &zvel, std::unique_ptr< MultiFab > &z_t_mf, const MultiFab &cc_src, const MultiFab &xmom_src, const MultiFab &ymom_src, const MultiFab &zmom_src, const MultiFab &buoyancy, const MultiFab *zmom_crse_rhs, Vector< std::unique_ptr< MultiFab >> &Tau_lev, MultiFab *SmnSmn, MultiFab *eddyDiffs, MultiFab *Hfx1, MultiFab *Hfx2, MultiFab *Hfx3, MultiFab *Q1fx1, MultiFab *Q1fx2, MultiFab *Q1fx3, MultiFab *Q2fx3, MultiFab *Diss, const Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, const Gpu::DeviceVector< BCRec > &domain_bcs_type_d, const Vector< BCRec > &domain_bcs_type_h, const MultiFab &z_phys_nd, const MultiFab &z_phys_cc, const MultiFab &ax, const MultiFab &ay, const MultiFab &az, const MultiFab &detJ, Gpu::DeviceVector< Real > &stretched_dz_d, Vector< MultiFab > &gradp, Vector< std::unique_ptr< MultiFab >> &mapfac, const eb_ &ebfact, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine)
Definition: ERF_SlowRhsPre.cpp:65
auto slow_rhs_fun_pre
Definition: ERF_TI_slow_rhs_pre.H:6
auto cons_to_prim
Definition: ERF_TI_utils.H:4
auto make_qt
Definition: ERF_TI_utils.H:52
auto update_terrain_stage
Definition: ERF_TI_utils.H:164
@ num_comps
Definition: ERF_IndexDefines.H:68
@ pi0_comp
Definition: ERF_IndexDefines.H:65
@ p0_comp
Definition: ERF_IndexDefines.H:64
@ th0_comp
Definition: ERF_IndexDefines.H:66
@ r0_comp
Definition: ERF_IndexDefines.H:63
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ qt
Definition: ERF_Kessler.H:27
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40