ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 int anelastic);
26 
27 void make_gradp_pert (int level,
28  const SolverChoice& solverChoice,
29  const amrex::Geometry& geom,
30  amrex::Vector<amrex::MultiFab>& S_data,
31  amrex::MultiFab& p0,
32  amrex::MultiFab& z_phys_nd,
33  amrex::MultiFab& z_phys_cc,
34  amrex::BCRec const* bcrec_ptr,
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  amrex::BCRec const* bcrec_ptr,
43  const eb_& ebfact,
44  amrex::Vector<amrex::MultiFab>& gradp,
45  const SolverChoice& solverChoice);
46 
47 void make_sources (int level, int nrk,
48  amrex::Real dt, amrex::Real time,
49  const amrex::Vector<amrex::MultiFab>& S_data,
50  const amrex::MultiFab& S_prim,
51  amrex::MultiFab& cc_source,
52  const amrex::MultiFab& base_state,
53  const amrex::MultiFab* z_phys_cc,
54  const amrex::MultiFab* qheating_rates,
55  amrex::MultiFab* terrain_blank,
56  const amrex::Geometry geom,
57  const SolverChoice& solverChoice,
58  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& mapfac,
59  const amrex::Real* dptr_rhotheta_src,
60  const amrex::Real* dptr_rhoqt_src,
61  const amrex::Real* dptr_wbar_sub,
62  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
63  InputSoundingData& input_sounding_data,
64  TurbulentPerturbation& turbPert,
65  bool is_slow_step);
66 
67 void make_mom_sources (amrex::Real time,
68  const amrex::Vector<amrex::MultiFab>& S_data,
69  amrex::MultiFab& z_phys_nd,
70  amrex::MultiFab& z_phys_cc,
71  const amrex::MultiFab& xvel,
72  const amrex::MultiFab& yvel,
73  const amrex::MultiFab& zvel,
74  amrex::MultiFab& xmom_source,
75  amrex::MultiFab& ymom_source,
76  amrex::MultiFab& zmom_source,
77  const amrex::MultiFab& base_state,
78  amrex::MultiFab* forest_drag,
79  amrex::MultiFab* terrain_blank,
80  amrex::MultiFab* cosPhi_mf,
81  amrex::MultiFab* sinPhi_mf,
82  const amrex::Geometry geom,
83  const SolverChoice& solverChoice,
84  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& mapfac,
85  const amrex::Real* dptr_rhotheta_src,
86  const amrex::Real* dptr_rhoqt_src,
87  const amrex::Real* dptr_wbar_sub,
88  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
89  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev,
90  InputSoundingData& input_sounding_data,
91  bool is_slow_step);
92 
93 void add_thin_body_sources (amrex::MultiFab& xmom_source,
94  amrex::MultiFab& ymom_source,
95  amrex::MultiFab& zmom_source,
96  std::unique_ptr<amrex::iMultiFab>& xflux_imask_lev,
97  std::unique_ptr<amrex::iMultiFab>& yflux_imask_lev,
98  std::unique_ptr<amrex::iMultiFab>& zflux_imask_lev,
99  std::unique_ptr<amrex::MultiFab>& thin_xforce_lev,
100  std::unique_ptr<amrex::MultiFab>& thin_yforce_lev,
101  std::unique_ptr<amrex::MultiFab>& thin_zforce_lev);
102 
103 #if defined(ERF_USE_NETCDF)
104 void
105 moist_set_rhs (const amrex::Box& tbx,
106  const amrex::Array4<amrex::Real const>& old_cons,
107  const amrex::Array4<amrex::Real const>& new_cons,
108  const amrex::Array4<amrex::Real >& cell_rhs,
109  const amrex::Real& bdy_time_interval,
110  const amrex::Real& start_bdy_time,
111  const amrex::Real& new_stage_time,
112  const amrex::Real& dt,
113  int width, int set_width,
114  const amrex::Box& domain,
115  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
116  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
117  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
118  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
119 #endif
120 
121 void ApplySpongeZoneBCsForCC (const SpongeChoice& spongeChoice,
122  const amrex::Geometry geom,
123  const amrex::Box& bx,
124  const amrex::Array4<amrex::Real>& cell_rhs,
125  const amrex::Array4<const amrex::Real>& cell_data,
126  const amrex::Array4<const amrex::Real>& r0,
127  const amrex::Array4<const amrex::Real>& z_phys_cc);
128 
129 void ApplySpongeZoneBCsForMom (const SpongeChoice& spongeChoice,
130  const amrex::Geometry geom,
131  const amrex::Box& tbx,
132  const amrex::Box& tby,
133  const amrex::Box& tbz,
134  const amrex::Array4<amrex::Real>& rho_u_rhs,
135  const amrex::Array4<amrex::Real>& rho_v_rhs,
136  const amrex::Array4<amrex::Real>& rho_w_rhs,
137  const amrex::Array4<const amrex::Real>& rho_u,
138  const amrex::Array4<const amrex::Real>& rho_v,
139  const amrex::Array4<const amrex::Real>& rho_w,
140  const amrex::Array4<const amrex::Real>& r0,
141  const amrex::Array4<const amrex::Real>& z_phys_nd,
142  const amrex::Array4<const amrex::Real>& z_phys_cc);
143 
145  const amrex::Geometry geom,
146  const amrex::Box& tbx,
147  const amrex::Box& tby,
148  const amrex::Array4<const amrex::Real>& cell_data,
149  const amrex::Array4<const amrex::Real>& z_phys_cc,
150  const amrex::Array4<amrex::Real>& rho_u_rhs,
151  const amrex::Array4<amrex::Real>& rho_v_rhs,
152  const amrex::Array4<const amrex::Real>& rho_u,
153  const amrex::Array4<const amrex::Real>& rho_v,
154  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev);
155 #endif
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, amrex::BCRec const *bcrec_ptr, const eb_ &ebfact, amrex::Vector< amrex::MultiFab > &gradp)
void compute_gradp(const amrex::MultiFab &p, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, amrex::BCRec const *bcrec_ptr, const eb_ &ebfact, amrex::Vector< amrex::MultiFab > &gradp, const SolverChoice &solverChoice)
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, const amrex::Vector< amrex::MultiFab > &S_data, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc, 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, InputSoundingData &input_sounding_data, bool is_slow_step)
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 int anelastic)
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:99
Definition: ERF_SpongeStruct.H:15
Definition: ERF_TurbPertStruct.H:22