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 "ERF_DataStruct.H"
7 #include "ERF_TurbPertStruct.H"
8 
9 #ifdef ERF_USE_EB
10 #include <AMReX_EBMultiFabUtil.H>
11 #endif
12 
13 /**
14  * Function for computing the buoyancy term to be used in the evolution
15  * equation for the z-component of momentum in the slow integrator. There
16  * are three options for how buoyancy is computed (two are the same in the absence of moisture).
17  */
18 void make_buoyancy (amrex::Vector< amrex::MultiFab>& S_data,
19  const amrex::MultiFab & S_prim,
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_sources (int level, int nrk,
28  amrex::Real dt, amrex::Real time,
29  amrex::Vector<amrex::MultiFab>& S_data,
30  const amrex::MultiFab& S_prim,
31  amrex::MultiFab& cc_source,
32  std::unique_ptr<amrex::MultiFab>& z_phys_cc,
33 #if defined(ERF_USE_RRTMGP)
34  const amrex::MultiFab* qheating_rates,
35 #endif
36  const amrex::Geometry geom,
37  const SolverChoice& solverChoice,
38  std::unique_ptr<amrex::MultiFab>& mapfac_u,
39  std::unique_ptr<amrex::MultiFab>& mapfac_v,
40  std::unique_ptr<amrex::MultiFab>& mapfac_m,
41  const amrex::Real* dptr_rhotheta_src,
42  const amrex::Real* dptr_rhoqt_src,
43  const amrex::Real* dptr_wbar_sub,
44  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
45  InputSoundingData& input_sounding_data,
46  TurbulentPerturbation& turbPert);
47 
48 void make_mom_sources (int level, int nrk,
49  amrex::Real dt, amrex::Real time,
50  amrex::Vector<amrex::MultiFab>& S_data,
51  const amrex::MultiFab& S_prim,
52  std::unique_ptr<amrex::MultiFab>& z_phys_nd,
53  std::unique_ptr<amrex::MultiFab>& z_phys_cc,
54  const amrex::MultiFab& xvel,
55  const amrex::MultiFab& yvel,
56  const amrex::MultiFab& wvel,
57  amrex::MultiFab& xmom_source,
58  amrex::MultiFab& ymom_source,
59  amrex::MultiFab& zmom_source,
60  const amrex::MultiFab& base_state,
61  amrex::MultiFab* forest_drag,
62  amrex::MultiFab* terrain_blank,
63  const amrex::Geometry geom,
64  const SolverChoice& solverChoice,
65  std::unique_ptr<amrex::MultiFab>& mapfac_m,
66  std::unique_ptr<amrex::MultiFab>& mapfac_u,
67  std::unique_ptr<amrex::MultiFab>& mapfac_v,
68  const amrex::Real* dptr_rhotheta_src,
69  const amrex::Real* dptr_rhoqt_src,
70  const amrex::Real* dptr_wbar_sub,
71  const amrex::Vector<amrex::Real*> d_rayleigh_ptrs_at_lev,
72  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev,
73  InputSoundingData& input_sounding_data,
74  const int n_qstate);
75 
76 void add_thin_body_sources (amrex::MultiFab& xmom_source,
77  amrex::MultiFab& ymom_source,
78  amrex::MultiFab& zmom_source,
79  std::unique_ptr<amrex::iMultiFab>& xflux_imask_lev,
80  std::unique_ptr<amrex::iMultiFab>& yflux_imask_lev,
81  std::unique_ptr<amrex::iMultiFab>& zflux_imask_lev,
82  std::unique_ptr<amrex::MultiFab>& thin_xforce_lev,
83  std::unique_ptr<amrex::MultiFab>& thin_yforce_lev,
84  std::unique_ptr<amrex::MultiFab>& thin_zforce_lev);
85 
86 #if defined(ERF_USE_NETCDF)
87 void
88 moist_set_rhs (const amrex::Box& tbx,
89  const amrex::Box& gtbx,
90  const amrex::Array4<amrex::Real const>& old_cons,
91  const amrex::Array4<amrex::Real const>& new_cons,
92  const amrex::Array4<amrex::Real >& cell_rhs,
93  const amrex::Real& bdy_time_interval,
94  const amrex::Real& start_bdy_time,
95  const amrex::Real& new_stage_time,
96  const amrex::Real& dt,
97  int width, int set_width,
98  const amrex::Box& domain,
99  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xlo,
100  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_xhi,
101  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_ylo,
102  amrex::Vector<amrex::Vector<amrex::FArrayBox>>& bdy_data_yhi);
103 #endif
104 
105 void ApplySpongeZoneBCsForCC (const SpongeChoice& spongeChoice,
106  const amrex::Geometry geom,
107  const amrex::Box& bx,
108  const amrex::Array4<amrex::Real>& cell_rhs,
109  const amrex::Array4<const amrex::Real>& cell_data);
110 
111 void ApplySpongeZoneBCsForMom (const SpongeChoice& spongeChoice,
112  const amrex::Geometry geom,
113  const amrex::Box& tbx,
114  const amrex::Box& tby,
115  const amrex::Box& tbz,
116  const amrex::Array4<amrex::Real>& rho_u_rhs,
117  const amrex::Array4<amrex::Real>& rho_v_rhs,
118  const amrex::Array4<amrex::Real>& rho_w_rhs,
119  const amrex::Array4<const amrex::Real>& rho_u,
120  const amrex::Array4<const amrex::Real>& rho_v,
121  const amrex::Array4<const amrex::Real>& rho_w);
122 
124  const amrex::Geometry geom,
125  const amrex::Box& tbx,
126  const amrex::Box& tby,
127  const amrex::Array4<const amrex::Real>& cell_data,
128  const amrex::Array4<amrex::Real>& rho_u_rhs,
129  const amrex::Array4<amrex::Real>& rho_v_rhs,
130  const amrex::Array4<const amrex::Real>& rho_u,
131  const amrex::Array4<const amrex::Real>& rho_v,
132  const amrex::Vector<amrex::Real*> d_sponge_ptrs_at_lev);
133 #endif
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)
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)
void make_buoyancy(amrex::Vector< amrex::MultiFab > &S_data, const amrex::MultiFab &S_prim, amrex::MultiFab &buoyancy, const amrex::Geometry geom, const SolverChoice &solverChoice, const amrex::MultiFab &base_state, const int n_qstate, const int anelastic)
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(int level, int nrk, amrex::Real dt, amrex::Real time, amrex::Vector< amrex::MultiFab > &S_data, const amrex::MultiFab &S_prim, std::unique_ptr< amrex::MultiFab > &z_phys_nd, std::unique_ptr< amrex::MultiFab > &z_phys_cc, const amrex::MultiFab &xvel, const amrex::MultiFab &yvel, const amrex::MultiFab &wvel, 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, const amrex::Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< amrex::MultiFab > &mapfac_m, std::unique_ptr< amrex::MultiFab > &mapfac_u, std::unique_ptr< amrex::MultiFab > &mapfac_v, 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, const int n_qstate)
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< 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 make_sources(int level, int nrk, amrex::Real dt, amrex::Real time, amrex::Vector< amrex::MultiFab > &S_data, const amrex::MultiFab &S_prim, amrex::MultiFab &cc_source, std::unique_ptr< amrex::MultiFab > &z_phys_cc, const amrex::Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< amrex::MultiFab > &mapfac_u, std::unique_ptr< amrex::MultiFab > &mapfac_v, std::unique_ptr< amrex::MultiFab > &mapfac_m, 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)
@ xvel
Definition: ERF_IndexDefines.H:130
@ yvel
Definition: ERF_IndexDefines.H:131
Definition: ERF_InputSoundingData.H:22
Definition: ERF_DataStruct.H:82
Definition: ERF_SpongeStruct.H:15
Definition: ERF_TurbPertStruct.H:18