ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SrcHeaders.H
Go to the documentation of this file.
1 #ifndef ERF_INTEGRATION_H_
2 #define ERF_INTEGRATION_H_
3 
4 #include <AMReX_MultiFab.H>
5 #include <AMReX_EBMultiFabUtil.H>
6 
7 #include "ERF_DataStruct.H"
9 #include "ERF_TurbPertStruct.H"
10 #include "ERF_EB.H"
11 
12 /**
13  * Function for computing the buoyancy term to be used in the evolution
14  * equation for the z-component of momentum in the slow integrator. There
15  * are three options for how buoyancy is computed (two are the same in the absence of moisture).
16  */
17 void make_buoyancy (const amrex::Vector<amrex::MultiFab>& S_data,
18  const amrex::MultiFab & S_prim,
19  const amrex::MultiFab & qt,
20  amrex::MultiFab& buoyancy,
21  const amrex::Geometry geom,
22  const SolverChoice& solverChoice,
23  const amrex::MultiFab& base_state,
24  const int n_qstate,
25  const eb_& ebfact,
26  const int anelastic);
27 
28 void make_gradp_pert (int level,
29  const SolverChoice& solverChoice,
30  const amrex::Geometry& geom,
31  amrex::Vector<amrex::MultiFab>& S_data,
32  amrex::MultiFab& p0,
33  amrex::MultiFab& z_phys_nd,
34  amrex::MultiFab& z_phys_cc,
35  const eb_& ebfact,
36  amrex::Vector<amrex::MultiFab>& gradp);
37 
38 void compute_gradp (const amrex::MultiFab& p,
39  const amrex::Geometry& geom,
40  amrex::MultiFab& z_phys_nd,
41  amrex::MultiFab& z_phys_cc,
42  const eb_& ebfact,
43  amrex::Vector<amrex::MultiFab>& gradp,
44  const SolverChoice& solverChoice);
45 
46 void compute_gradp_interpz (const amrex::MultiFab& p,
47  const amrex::Geometry& geom,
48  amrex::MultiFab& z_phys_nd,
49  amrex::MultiFab& z_phys_cc,
50  amrex::Vector<amrex::MultiFab>& gradp,
51  const SolverChoice& solverChoice);
52 
53 void make_sources (int level, int nrk,
54  amrex::Real dt, amrex::Real time,
55  const amrex::Vector<amrex::MultiFab>& S_data,
56  const amrex::MultiFab& S_prim,
57  amrex::MultiFab& cc_source,
58  const amrex::MultiFab& base_state,
59  const amrex::MultiFab* z_phys_cc,
60  const amrex::MultiFab* qheating_rates,
61  amrex::MultiFab* terrain_blank,
62  const amrex::Geometry geom,
63  const SolverChoice& solverChoice,
64  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& mapfac,
65  const amrex::Real* dptr_rhotheta_src,
66  const amrex::Real* dptr_rhoqt_src,
67  const amrex::Real* dptr_wbar_sub,
68  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
69  InputSoundingData& input_sounding_data,
70  TurbulentPerturbation& turbPert,
71  bool is_slow_step);
72 
74  amrex::Real dt,
75  const amrex::Vector<amrex::MultiFab>& S_data,
76  amrex::MultiFab& z_phys_nd,
77  amrex::MultiFab& z_phys_cc,
78  amrex::Vector<amrex::Real>& stretched_dz_h,
79  const amrex::MultiFab& xvel,
80  const amrex::MultiFab& yvel,
81  const amrex::MultiFab& zvel,
82  amrex::MultiFab& xmom_source,
83  amrex::MultiFab& ymom_source,
84  amrex::MultiFab& zmom_source,
85  const amrex::MultiFab& base_state,
86  amrex::MultiFab* forest_drag,
87  amrex::MultiFab* terrain_blank,
88  amrex::MultiFab* cosPhi_mf,
89  amrex::MultiFab* sinPhi_mf,
90  const amrex::Geometry geom,
91  const SolverChoice& solverChoice,
92  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& mapfac,
93  const amrex::Real* dptr_rhotheta_src,
94  const amrex::Real* dptr_rhoqt_src,
95  const amrex::Real* dptr_wbar_sub,
96  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
97  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev,
98  const amrex::Vector<amrex::MultiFab>* forecast_state_at_lev,
99  InputSoundingData& input_sounding_data,
100  bool is_slow_step);
101 
102 void add_thin_body_sources (amrex::MultiFab& xmom_source,
103  amrex::MultiFab& ymom_source,
104  amrex::MultiFab& zmom_source,
105  std::unique_ptr<amrex::iMultiFab>& xflux_imask_lev,
106  std::unique_ptr<amrex::iMultiFab>& yflux_imask_lev,
107  std::unique_ptr<amrex::iMultiFab>& zflux_imask_lev,
108  std::unique_ptr<amrex::MultiFab>& thin_xforce_lev,
109  std::unique_ptr<amrex::MultiFab>& thin_yforce_lev,
110  std::unique_ptr<amrex::MultiFab>& thin_zforce_lev);
111 
112 #if defined(ERF_USE_NETCDF)
113 void
114 moist_set_rhs (const amrex::Geometry& geom,
115  const amrex::Box& tbx,
116  const amrex::Array4<amrex::Real const>& old_cons,
117  const amrex::Array4<amrex::Real const>& new_cons,
118  const amrex::Array4<amrex::Real >& cell_rhs,
119  const amrex::Real& bdy_time_interval,
120  const amrex::Real& new_stage_time,
121  const amrex::Real& dt,
122  const amrex::Real& stop_time_elapsed,
123  int width, int set_width,
124  const amrex::Box& domain,
125  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
126  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
127  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
128  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
129 #endif
130 
131 void ApplySpongeZoneBCsForCC (const SpongeChoice& spongeChoice,
132  const amrex::Geometry geom,
133  const amrex::Box& bx,
134  const amrex::Array4<amrex::Real>& cell_rhs,
135  const amrex::Array4<const amrex::Real>& cell_data,
136  const amrex::Array4<const amrex::Real>& r0,
137  const amrex::Array4<const amrex::Real>& z_phys_cc);
138 
139 void ApplySpongeZoneBCsForMom (const SpongeChoice& spongeChoice,
140  const amrex::Geometry geom,
141  const amrex::Box& tbx,
142  const amrex::Box& tby,
143  const amrex::Box& tbz,
144  const amrex::Array4<amrex::Real>& rho_u_rhs,
145  const amrex::Array4<amrex::Real>& rho_v_rhs,
146  const amrex::Array4<amrex::Real>& rho_w_rhs,
147  const amrex::Array4<const amrex::Real>& rho_u,
148  const amrex::Array4<const amrex::Real>& rho_v,
149  const amrex::Array4<const amrex::Real>& rho_w,
150  const amrex::Array4<const amrex::Real>& r0,
151  const amrex::Array4<const amrex::Real>& z_phys_nd,
152  const amrex::Array4<const amrex::Real>& z_phys_cc);
153 
155  const amrex::Geometry geom,
156  const amrex::Box& tbx,
157  const amrex::Box& tby,
158  const amrex::Array4<const amrex::Real>& cell_data,
159  const amrex::Array4<const amrex::Real>& z_phys_cc,
160  const amrex::Array4<amrex::Real>& rho_u_rhs,
161  const amrex::Array4<amrex::Real>& rho_v_rhs,
162  const amrex::Array4<const amrex::Real>& rho_u,
163  const amrex::Array4<const amrex::Real>& rho_v,
164  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev);
165 
166 void ApplyBndryForcing_Forecast (const SolverChoice& solverChoice,
167  const amrex::Geometry geom,
168  const amrex::Box& tbx,
169  const amrex::Box& tby,
170  const amrex::Box& tbz,
171  const amrex::Array4<const amrex::Real>& z_phys_nd,
172  const amrex::Array4<amrex::Real>& rho_u_rhs,
173  const amrex::Array4<amrex::Real>& rho_v_rhs,
174  const amrex::Array4<amrex::Real>& rho_w_rhs,
175  const amrex::Array4<const amrex::Real>& rho_u,
176  const amrex::Array4<const amrex::Real>& rho_v,
177  const amrex::Array4<const amrex::Real>& rho_w,
178  const amrex::Array4<const amrex::Real>& rho_u_initial_state,
179  const amrex::Array4<const amrex::Real>& rho_v_initial_state,
180  const amrex::Array4<const amrex::Real>& rho_w_initial_state,
181  const amrex::Array4<const amrex::Real>& cons_initial_state);
182 
183 #endif
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void make_gradp_pert(int level, const SolverChoice &solverChoice, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_data, amrex::MultiFab &p0, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, const eb_ &ebfact, amrex::Vector< amrex::MultiFab > &gradp)
void ApplyBndryForcing_Forecast(const SolverChoice &solverChoice, const amrex::Geometry geom, const amrex::Box &tbx, const amrex::Box &tby, const amrex::Box &tbz, const amrex::Array4< const amrex::Real > &z_phys_nd, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &rho_w, const amrex::Array4< const amrex::Real > &rho_u_initial_state, const amrex::Array4< const amrex::Real > &rho_v_initial_state, const amrex::Array4< const amrex::Real > &rho_w_initial_state, const amrex::Array4< const amrex::Real > &cons_initial_state)
void add_thin_body_sources(amrex::MultiFab &xmom_source, amrex::MultiFab &ymom_source, amrex::MultiFab &zmom_source, std::unique_ptr< amrex::iMultiFab > &xflux_imask_lev, std::unique_ptr< amrex::iMultiFab > &yflux_imask_lev, std::unique_ptr< amrex::iMultiFab > &zflux_imask_lev, std::unique_ptr< amrex::MultiFab > &thin_xforce_lev, std::unique_ptr< amrex::MultiFab > &thin_yforce_lev, std::unique_ptr< amrex::MultiFab > &thin_zforce_lev)
void make_mom_sources(amrex::Real time, amrex::Real dt, const amrex::Vector< amrex::MultiFab > &S_data, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, amrex::Vector< amrex::Real > &stretched_dz_h, const amrex::MultiFab &xvel, const amrex::MultiFab &yvel, const amrex::MultiFab &zvel, amrex::MultiFab &xmom_source, amrex::MultiFab &ymom_source, amrex::MultiFab &zmom_source, const amrex::MultiFab &base_state, amrex::MultiFab *forest_drag, amrex::MultiFab *terrain_blank, amrex::MultiFab *cosPhi_mf, amrex::MultiFab *sinPhi_mf, const amrex::Geometry geom, const SolverChoice &solverChoice, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &mapfac, const amrex::Real *dptr_rhotheta_src, const amrex::Real *dptr_rhoqt_src, const amrex::Real *dptr_wbar_sub, const amrex::Vector< amrex::Real * > d_rayleigh_ptrs_at_lev, const amrex::Vector< amrex::Real * > d_sponge_ptrs_at_lev, const amrex::Vector< amrex::MultiFab > *forecast_state_at_lev, InputSoundingData &input_sounding_data, bool is_slow_step)
void compute_gradp_interpz(const amrex::MultiFab &p, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, amrex::Vector< amrex::MultiFab > &gradp, const SolverChoice &solverChoice)
void make_buoyancy(const amrex::Vector< amrex::MultiFab > &S_data, const amrex::MultiFab &S_prim, const amrex::MultiFab &qt, amrex::MultiFab &buoyancy, const amrex::Geometry geom, const SolverChoice &solverChoice, const amrex::MultiFab &base_state, const int n_qstate, const eb_ &ebfact, const int anelastic)
void compute_gradp(const amrex::MultiFab &p, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, const eb_ &ebfact, amrex::Vector< amrex::MultiFab > &gradp, const SolverChoice &solverChoice)
void ApplySpongeZoneBCsForMom_ReadFromFile(const SpongeChoice &spongeChoice, const amrex::Geometry geom, const amrex::Box &tbx, const amrex::Box &tby, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &z_phys_cc, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Vector< amrex::Real * > d_sponge_ptrs_at_lev)
void ApplySpongeZoneBCsForCC(const SpongeChoice &spongeChoice, const amrex::Geometry geom, const amrex::Box &bx, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &r0, const amrex::Array4< const amrex::Real > &z_phys_cc)
void ApplySpongeZoneBCsForMom(const SpongeChoice &spongeChoice, const amrex::Geometry geom, const amrex::Box &tbx, const amrex::Box &tby, const amrex::Box &tbz, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &rho_w, const amrex::Array4< const amrex::Real > &r0, const amrex::Array4< const amrex::Real > &z_phys_nd, const amrex::Array4< const amrex::Real > &z_phys_cc)
void make_sources(int level, int nrk, amrex::Real dt, amrex::Real time, const amrex::Vector< amrex::MultiFab > &S_data, const amrex::MultiFab &S_prim, amrex::MultiFab &cc_source, const amrex::MultiFab &base_state, const amrex::MultiFab *z_phys_cc, const amrex::MultiFab *qheating_rates, amrex::MultiFab *terrain_blank, const amrex::Geometry geom, const SolverChoice &solverChoice, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &mapfac, const amrex::Real *dptr_rhotheta_src, const amrex::Real *dptr_rhoqt_src, const amrex::Real *dptr_wbar_sub, const amrex::Vector< amrex::Real * > d_rayleigh_ptrs_at_lev, InputSoundingData &input_sounding_data, TurbulentPerturbation &turbPert, bool is_slow_step)
Definition: ERF_EB.H:13
@ qt
Definition: ERF_Kessler.H:27
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
real(c_double), parameter p0
Definition: ERF_module_model_constants.F90:40
Definition: ERF_InputSoundingData.H:22
Definition: ERF_DataStruct.H:123
Definition: ERF_SpongeStruct.H:15
Definition: ERF_TurbPertStruct.H:22