ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_BuoyancyUtils.H
Go to the documentation of this file.
1 #ifndef ERF_BUOYANCY_UTILS_H_
2 #define ERF_BUOYANCY_UTILS_H_
3 
4 #include <ERF_EOS.H>
5 #include <ERF_Constants.H>
6 
7 AMREX_GPU_DEVICE
8 AMREX_FORCE_INLINE
9 amrex::Real
10 buoyancy_dry_anelastic (int& i, int& j, int& k,
11  amrex::Real const& grav_gpu,
12  const amrex::Array4<const amrex::Real>& r0_arr,
13  const amrex::Array4<const amrex::Real>& th0_arr,
14  const amrex::Array4<const amrex::Real>& cell_data)
15 {
16  // Note: this is the same term as the moist anelastic buoyancy when qv = qc = qt = 0
17  amrex::Real theta_d_hi = cell_data(i,j,k ,RhoTheta_comp)/cell_data(i,j,k ,Rho_comp);
18  amrex::Real theta_d_lo = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
19 
20  amrex::Real theta_d_wface = amrex::Real(0.5) * (theta_d_lo + theta_d_hi);
21  amrex::Real theta_d0_wface = amrex::Real(0.5) * (th0_arr(i,j,k) + th0_arr(i,j,k-1));
22  amrex::Real rho0_wface = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
23 
24  return (-rho0_wface * grav_gpu * (theta_d_wface - theta_d0_wface) / theta_d0_wface);
25 }
26 
27 AMREX_GPU_DEVICE
28 AMREX_FORCE_INLINE
29 amrex::Real
30 buoyancy_moist_anelastic (int& i, int& j, int& k,
31  amrex::Real const& grav_gpu,
32  amrex::Real const& rv_over_rd,
33  const amrex::Array4<const amrex::Real>& r0_arr,
34  const amrex::Array4<const amrex::Real>& th0_arr,
35  const amrex::Array4<const amrex::Real>& qv0_arr,
36  const amrex::Array4<const amrex::Real>& cell_data,
37  const amrex::Array4<const amrex::Real>& qt_arr)
38 {
39  amrex::Real Fact = (1.0 - rv_over_rd); // = -0.61
40 
41  amrex::Real th_d_hi = cell_data(i,j,k ,RhoTheta_comp)/cell_data(i,j,k ,Rho_comp);
42  amrex::Real th_d_lo = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
43 
44  amrex::Real qv_hi = cell_data(i,j,k ,RhoQ1_comp) /cell_data(i,j,k ,Rho_comp);
45  amrex::Real qv_lo = cell_data(i,j,k-1,RhoQ1_comp) /cell_data(i,j,k-1,Rho_comp);
46 
47  amrex::Real qt_hi = qt_arr(i,j,k ) - qv_hi;
48  amrex::Real th_v_hi = th_d_hi * (1.0 - Fact*qt_hi - rv_over_rd*qt_hi);
49 
50  amrex::Real qt_lo = qt_arr(i,j,k-1) - qv_lo;
51  amrex::Real th_v_lo = th_d_lo * (1.0 - Fact*qt_lo - rv_over_rd*qt_lo);
52 
53  amrex::Real th_v0_hi = th0_arr(i,j,k ) * (1.0 - Fact*qv0_arr(i,j,k ));
54  amrex::Real th_v0_lo = th0_arr(i,j,k-1) * (1.0 - Fact*qv0_arr(i,j,k-1));
55 
56  amrex::Real qv0_hi = qv0_arr(i,j,k );
57  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
58 
59  amrex::Real q_hi = (qv_hi-qv0_hi) - qt_hi + (th_v_hi - th_v0_hi) / th_v0_hi;
60  amrex::Real q_lo = (qv_lo-qv0_lo) - qt_lo + (th_v_lo - th_v0_lo) / th_v0_lo;
61 
62  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
63 
64  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
65 
66  return (-r0avg * grav_gpu * qavg);
67 }
68 
69 AMREX_GPU_DEVICE
70 AMREX_FORCE_INLINE
71 amrex::Real
72 buoyancy_rhopert (int& i, int& j, int& k,
73  amrex::Real const& grav_gpu,
74  const amrex::Array4<const amrex::Real>& r0_arr,
75  const amrex::Array4<const amrex::Real>& cell_data,
76  const amrex::Array4<const amrex::Real>& qt_arr)
77 {
78  amrex::Real rhop_hi = cell_data(i,j,k ,Rho_comp) * (1.0 + qt_arr(i,j,k ))- r0_arr(i,j,k );
79  amrex::Real rhop_lo = cell_data(i,j,k-1,Rho_comp) * (1.0 + qt_arr(i,j,k-1))- r0_arr(i,j,k-1);
80  return( grav_gpu * amrex::Real(0.5) * ( rhop_hi + rhop_lo ) );
81 }
82 
83 AMREX_GPU_DEVICE
84 AMREX_FORCE_INLINE
85 amrex::Real
86 buoyancy_dry_Tpert (int& i, int& j, int& k,
87  amrex::Real const& grav_gpu,
88  amrex::Real const& rd_over_cp,
89  const amrex::Array4<const amrex::Real>& r0_arr,
90  const amrex::Array4<const amrex::Real>& p0_arr,
91  const amrex::Array4<const amrex::Real>& th0_arr,
92  const amrex::Array4<const amrex::Real>& cell_data)
93 {
94  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k), th0_arr(i,j,k), rd_over_cp);
95  amrex::Real t_hi = getTgivenRandRTh(cell_data(i,j,k ,Rho_comp), cell_data(i,j,k ,RhoTheta_comp));
96 
97  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), th0_arr(i,j,k-1), rd_over_cp);
98  amrex::Real t_lo = getTgivenRandRTh(cell_data(i,j,k-1,Rho_comp), cell_data(i,j,k-1,RhoTheta_comp));
99 
100  amrex::Real tprime_hi = (t_hi-t0_hi)/t0_hi;
101  amrex::Real tprime_lo = (t_lo-t0_lo)/t0_lo;
102 
103  amrex::Real tp_avg = amrex::Real(0.5) * (tprime_hi + tprime_lo);
104  amrex::Real r0_avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
105 
106  return ( -r0_avg * grav_gpu * tp_avg);
107 }
108 
109 AMREX_GPU_DEVICE
110 AMREX_FORCE_INLINE
111 amrex::Real
112 buoyancy_dry_Thpert (int& i, int& j, int& k,
113  amrex::Real const& grav_gpu,
114  const amrex::Array4<const amrex::Real>& r0_arr,
115  const amrex::Array4<const amrex::Real>& th0_arr,
116  const amrex::Array4<const amrex::Real>& cell_prim)
117 {
118  amrex::Real thetaprime_hi = (cell_prim(i,j,k ,PrimTheta_comp) - th0_arr(i,j,k )) / th0_arr(i,j,k );
119  amrex::Real thetaprime_lo = (cell_prim(i,j,k-1,PrimTheta_comp) - th0_arr(i,j,k-1)) / th0_arr(i,j,k-1);
120 
121  amrex::Real thp_avg = amrex::Real(0.5) * (thetaprime_hi + thetaprime_lo);
122  amrex::Real r0_avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
123 
124  return ( -r0_avg * grav_gpu * thp_avg);
125 }
126 
127 AMREX_GPU_DEVICE
128 AMREX_FORCE_INLINE
129 amrex::Real
130 buoyancy_moist_Tpert (int& i, int& j, int& k,
131  const int& n_qstate,
132  amrex::Real const& grav_gpu,
133  amrex::Real const& rd_over_cp,
134  const amrex::Array4<const amrex::Real>& r0_arr,
135  const amrex::Array4<const amrex::Real>& th0_arr,
136  const amrex::Array4<const amrex::Real>& qv0_arr,
137  const amrex::Array4<const amrex::Real>& p0_arr,
138  const amrex::Array4<const amrex::Real>& cell_prim,
139  const amrex::Array4<const amrex::Real>& cell_data,
140  const amrex::Array4<const amrex::Real>& qt_arr)
141 {
142  amrex::Real qv_hi = (n_qstate >= 1) ? cell_prim(i,j,k ,PrimQ1_comp) : 0.0;
143  amrex::Real qv_lo = (n_qstate >= 1) ? cell_prim(i,j,k-1,PrimQ1_comp) : 0.0;
144 
145  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k ), th0_arr(i,j,k ), rd_over_cp);
146  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), th0_arr(i,j,k-1), rd_over_cp);
147 
148  amrex::Real dt_hi = getTgivenRandRTh(cell_data(i,j,k ,Rho_comp), cell_data(i,j,k ,RhoTheta_comp), qv_hi) - t0_hi;
149  amrex::Real dt_lo = getTgivenRandRTh(cell_data(i,j,k-1,Rho_comp), cell_data(i,j,k-1,RhoTheta_comp), qv_lo) - t0_lo;
150 
151  amrex::Real qv0_hi = qv0_arr(i,j,k );
152  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
153 
154  amrex::Real q_hi = amrex::Real(0.61) * (qv_hi-qv0_hi) - (qt_arr(i,j,k )-qv_hi) + dt_hi/t0_hi;
155  amrex::Real q_lo = amrex::Real(0.61) * (qv_lo-qv0_lo) - (qt_arr(i,j,k-1)-qv_lo) + dt_lo/t0_lo;
156 
157  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
158 
159  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
160 
161  return ( -r0avg * grav_gpu * qavg);
162 }
163 
164 AMREX_GPU_DEVICE
165 AMREX_FORCE_INLINE
166 amrex::Real
167 buoyancy_moist_Thpert (int& i, int& j, int& k,
168  const int& n_qstate,
169  amrex::Real const& grav_gpu,
170  const amrex::Array4<const amrex::Real>& r0_arr,
171  const amrex::Array4<const amrex::Real>& th0_arr,
172  const amrex::Array4<const amrex::Real>& qv0_arr,
173  const amrex::Array4<const amrex::Real>& cell_prim,
174  const amrex::Array4<const amrex::Real>& qt_arr)
175 {
176  amrex::Real qv_hi = (n_qstate >= 1) ? cell_prim(i,j,k ,PrimQ1_comp) : 0.0;
177  amrex::Real qv_lo = (n_qstate >= 1) ? cell_prim(i,j,k-1,PrimQ1_comp) : 0.0;
178 
179  amrex::Real qt_hi = qt_arr(i,j,k ) - qv_hi;
180  amrex::Real qt_lo = qt_arr(i,j,k-1) - qv_lo;
181 
182  amrex::Real dth_hi = cell_prim(i,j,k ,PrimTheta_comp) - th0_arr(i,j,k );
183  amrex::Real dth_lo = cell_prim(i,j,k-1,PrimTheta_comp) - th0_arr(i,j,k-1);
184 
185  amrex::Real qv0_hi = qv0_arr(i,j,k );
186  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
187 
188  amrex::Real q_hi = amrex::Real(0.61) * (qv_hi-qv0_hi) - qt_hi + dth_hi/th0_arr(i,j,k );
189  amrex::Real q_lo = amrex::Real(0.61) * (qv_lo-qv0_lo) - qt_lo + dth_lo/th0_arr(i,j,k-1);
190 
191  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
192 
193  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
194 
195  return ( -r0avg * grav_gpu * qavg);
196 }
197 
198 // **************************************************************************************
199 // Routines below this line are not currently used
200 // **************************************************************************************
201 
202 AMREX_GPU_DEVICE
203 AMREX_FORCE_INLINE
204 amrex::Real
205 buoyancy_dry_anelastic_T (int& i, int& j, int& k,
206  amrex::Real const& grav_gpu,
207  amrex::Real const& rd_over_cp,
208  const amrex::Array4<const amrex::Real>& r0_arr,
209  const amrex::Array4<const amrex::Real>& p0_arr,
210  const amrex::Array4<const amrex::Real>& cell_data)
211 {
212  amrex::Real rt0_hi = getRhoThetagivenP(p0_arr(i,j,k));
213  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k), rt0_hi/r0_arr(i,j,k), rd_over_cp);
214  amrex::Real t_hi = getTgivenPandTh(p0_arr(i,j,k), cell_data(i,j,k,RhoTheta_comp)/r0_arr(i,j,k), rd_over_cp);
215  amrex::Real q_hi = (t_hi-t0_hi)/t0_hi;
216 
217  amrex::Real rt0_lo = getRhoThetagivenP(p0_arr(i,j,k-1));
218  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), rt0_lo/r0_arr(i,j,k-1), rd_over_cp);
219  amrex::Real t_lo = getTgivenPandTh(p0_arr(i,j,k-1), cell_data(i,j,k-1,RhoTheta_comp)/r0_arr(i,j,k-1), rd_over_cp);
220  amrex::Real q_lo = (t_lo-t0_lo)/t0_lo;
221 
222  amrex::Real r0_q_avg = amrex::Real(0.5) * (r0_arr(i,j,k) * q_hi + r0_arr(i,j,k-1) * q_lo);
223  return (-r0_q_avg * grav_gpu);
224 }
225 
226 #endif
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_moist_anelastic(int &i, int &j, int &k, amrex::Real const &grav_gpu, amrex::Real const &rv_over_rd, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &qv0_arr, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &qt_arr)
Definition: ERF_BuoyancyUtils.H:30
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_moist_Tpert(int &i, int &j, int &k, const int &n_qstate, amrex::Real const &grav_gpu, amrex::Real const &rd_over_cp, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &qv0_arr, const amrex::Array4< const amrex::Real > &p0_arr, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &qt_arr)
Definition: ERF_BuoyancyUtils.H:130
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_moist_Thpert(int &i, int &j, int &k, const int &n_qstate, amrex::Real const &grav_gpu, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &qv0_arr, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &qt_arr)
Definition: ERF_BuoyancyUtils.H:167
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_dry_anelastic(int &i, int &j, int &k, amrex::Real const &grav_gpu, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &cell_data)
Definition: ERF_BuoyancyUtils.H:10
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_dry_Thpert(int &i, int &j, int &k, amrex::Real const &grav_gpu, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &cell_prim)
Definition: ERF_BuoyancyUtils.H:112
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_dry_Tpert(int &i, int &j, int &k, amrex::Real const &grav_gpu, amrex::Real const &rd_over_cp, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &p0_arr, const amrex::Array4< const amrex::Real > &th0_arr, const amrex::Array4< const amrex::Real > &cell_data)
Definition: ERF_BuoyancyUtils.H:86
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_dry_anelastic_T(int &i, int &j, int &k, amrex::Real const &grav_gpu, amrex::Real const &rd_over_cp, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &p0_arr, const amrex::Array4< const amrex::Real > &cell_data)
Definition: ERF_BuoyancyUtils.H:205
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_rhopert(int &i, int &j, int &k, amrex::Real const &grav_gpu, const amrex::Array4< const amrex::Real > &r0_arr, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &qt_arr)
Definition: ERF_BuoyancyUtils.H:72
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhoThetagivenP(const amrex::Real p, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:175
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenPandTh(const amrex::Real P, const amrex::Real th, const amrex::Real rdOcp)
Definition: ERF_EOS.H:32
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:53
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42