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 {
38  amrex::Real Fact = (1.0 - rv_over_rd); // = -0.61
39 
40  amrex::Real theta_d_hi = cell_data(i,j,k,RhoTheta_comp)/cell_data(i,j,k,Rho_comp);
41  amrex::Real qv_hi = cell_data(i,j,k,RhoQ1_comp) /cell_data(i,j,k,Rho_comp);
42  amrex::Real qc_hi = cell_data(i,j,k,RhoQ2_comp) /cell_data(i,j,k,Rho_comp);
43  amrex::Real qt_hi = qv_hi + qc_hi;
44  amrex::Real theta_v_hi = theta_d_hi * (1.0 - Fact*qt_hi - rv_over_rd*qc_hi);
45 
46  amrex::Real theta_d_lo = cell_data(i,j,k-1,RhoTheta_comp)/cell_data(i,j,k-1,Rho_comp);
47  amrex::Real qv_lo = cell_data(i,j,k-1,RhoQ1_comp) /cell_data(i,j,k-1,Rho_comp);
48  amrex::Real qc_lo = cell_data(i,j,k-1,RhoQ2_comp) /cell_data(i,j,k-1,Rho_comp);
49  amrex::Real qt_lo = qv_lo + qc_lo;
50  amrex::Real theta_v_lo = theta_d_lo * (1.0 - Fact*qt_lo - rv_over_rd*qc_lo);
51 
52  amrex::Real theta_v0_hi = th0_arr(i,j,k ) * (1.0 - Fact*qv0_arr(i,j,k ));
53  amrex::Real theta_v0_lo = th0_arr(i,j,k-1) * (1.0 - Fact*qv0_arr(i,j,k-1));
54  amrex::Real qv0_hi = qv0_arr(i,j,k );
55  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
56 
57  amrex::Real q_hi = (qv_hi-qv0_hi) - qc_hi + (theta_v_hi - theta_v0_hi) / theta_v0_hi;
58  amrex::Real q_lo = (qv_lo-qv0_lo) - qc_lo + (theta_v_lo - theta_v0_lo) / theta_v0_lo;
59  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
60 
61  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
62 
63  return (-r0avg * grav_gpu * qavg);
64 }
65 
66 AMREX_GPU_DEVICE
67 AMREX_FORCE_INLINE
68 amrex::Real
69 buoyancy_rhopert (int& i, int& j, int& k,
70  const int& n_qstate,
71  amrex::Real const& grav_gpu,
72  const amrex::Array4<const amrex::Real>& r0_arr,
73  const amrex::Array4<const amrex::Real>& cell_data)
74 {
75  amrex::Real rhop_hi = cell_data(i,j,k ,Rho_comp) - r0_arr(i,j,k );
76  amrex::Real rhop_lo = cell_data(i,j,k-1,Rho_comp) - r0_arr(i,j,k-1);
77  for (int q_offset(0); q_offset<n_qstate; ++q_offset) {
78  rhop_hi += cell_data(i,j,k ,RhoQ1_comp+q_offset);
79  rhop_lo += cell_data(i,j,k-1,RhoQ1_comp+q_offset);
80  }
81  return( grav_gpu * amrex::Real(0.5) * ( rhop_hi + rhop_lo ) );
82 }
83 
84 AMREX_GPU_DEVICE
85 AMREX_FORCE_INLINE
86 amrex::Real
87 buoyancy_dry_Tpert (int& i, int& j, int& k,
88  amrex::Real const& grav_gpu,
89  amrex::Real const& rd_over_cp,
90  const amrex::Array4<const amrex::Real>& r0_arr,
91  const amrex::Array4<const amrex::Real>& p0_arr,
92  const amrex::Array4<const amrex::Real>& th0_arr,
93  const amrex::Array4<const amrex::Real>& cell_data)
94 {
95  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k), th0_arr(i,j,k), rd_over_cp);
96  amrex::Real t_hi = getTgivenRandRTh(cell_data(i,j,k ,Rho_comp), cell_data(i,j,k ,RhoTheta_comp));
97 
98  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), th0_arr(i,j,k-1), rd_over_cp);
99  amrex::Real t_lo = getTgivenRandRTh(cell_data(i,j,k-1,Rho_comp), cell_data(i,j,k-1,RhoTheta_comp));
100 
101  amrex::Real tprime_hi = (t_hi-t0_hi)/t0_hi;
102  amrex::Real tprime_lo = (t_lo-t0_lo)/t0_lo;
103 
104  amrex::Real tp_avg = amrex::Real(0.5) * (tprime_hi + tprime_lo);
105  amrex::Real r0_avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
106 
107  return ( -r0_avg * grav_gpu * tp_avg);
108 }
109 
110 AMREX_GPU_DEVICE
111 AMREX_FORCE_INLINE
112 amrex::Real
113 buoyancy_dry_Thpert (int& i, int& j, int& k,
114  amrex::Real const& grav_gpu,
115  const amrex::Array4<const amrex::Real>& r0_arr,
116  const amrex::Array4<const amrex::Real>& th0_arr,
117  const amrex::Array4<const amrex::Real>& cell_prim)
118 {
119  amrex::Real thetaprime_hi = (cell_prim(i,j,k ,PrimTheta_comp) - th0_arr(i,j,k )) / th0_arr(i,j,k );
120  amrex::Real thetaprime_lo = (cell_prim(i,j,k-1,PrimTheta_comp) - th0_arr(i,j,k-1)) / th0_arr(i,j,k-1);
121 
122  amrex::Real thp_avg = amrex::Real(0.5) * (thetaprime_hi + thetaprime_lo);
123  amrex::Real r0_avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
124 
125  return ( -r0_avg * grav_gpu * thp_avg);
126 }
127 
128 AMREX_GPU_DEVICE
129 AMREX_FORCE_INLINE
130 amrex::Real
131 buoyancy_moist_Tpert (int& i, int& j, int& k,
132  const int& n_qstate,
133  amrex::Real const& grav_gpu,
134  amrex::Real const& rd_over_cp,
135  const amrex::Array4<const amrex::Real>& r0_arr,
136  const amrex::Array4<const amrex::Real>& th0_arr,
137  const amrex::Array4<const amrex::Real>& qv0_arr,
138  const amrex::Array4<const amrex::Real>& p0_arr,
139  const amrex::Array4<const amrex::Real>& cell_prim,
140  const amrex::Array4<const amrex::Real>& cell_data)
141 {
142  amrex::Real qv_hi = (n_qstate >= 1) ? cell_prim(i,j,k ,PrimQ1_comp) : 0.0;
143  amrex::Real qc_hi = (n_qstate >= 2) ? cell_prim(i,j,k ,PrimQ2_comp) : 0.0;
144  amrex::Real qp_hi = (n_qstate >= 3) ? cell_prim(i,j,k ,PrimQ3_comp) : 0.0;
145 
146  amrex::Real qv_lo = (n_qstate >= 1) ? cell_prim(i,j,k-1,PrimQ1_comp) : 0.0;
147  amrex::Real qc_lo = (n_qstate >= 2) ? cell_prim(i,j,k-1,PrimQ2_comp) : 0.0;
148  amrex::Real qp_lo = (n_qstate >= 3) ? cell_prim(i,j,k-1,PrimQ3_comp) : 0.0;
149 
150  amrex::Real t_hi = getTgivenRandRTh(cell_data(i,j,k ,Rho_comp), cell_data(i,j,k ,RhoTheta_comp), qv_hi);
151  amrex::Real t_lo = getTgivenRandRTh(cell_data(i,j,k-1,Rho_comp), cell_data(i,j,k-1,RhoTheta_comp), qv_lo);
152 
153  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k ), th0_arr(i,j,k ), rd_over_cp);
154  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), th0_arr(i,j,k-1), rd_over_cp);
155 
156  amrex::Real qv0_hi = qv0_arr(i,j,k );
157  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
158 
159  amrex::Real q_hi = amrex::Real(0.61) * (qv_hi-qv0_hi) - (qc_hi + qp_hi) + (t_hi-t0_hi)/t0_hi;
160  amrex::Real q_lo = amrex::Real(0.61) * (qv_lo-qv0_lo) - (qc_lo + qp_lo) + (t_lo-t0_lo)/t0_lo;
161  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
162 
163  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
164 
165  return ( -r0avg * grav_gpu * qavg);
166 }
167 
168 AMREX_GPU_DEVICE
169 AMREX_FORCE_INLINE
170 amrex::Real
171 buoyancy_moist_Thpert (int& i, int& j, int& k,
172  const int& n_qstate,
173  amrex::Real const& grav_gpu,
174  const amrex::Array4<const amrex::Real>& r0_arr,
175  const amrex::Array4<const amrex::Real>& th0_arr,
176  const amrex::Array4<const amrex::Real>& qv0_arr,
177  const amrex::Array4<const amrex::Real>& cell_prim)
178 {
179  amrex::Real qv_hi = (n_qstate >= 1) ? cell_prim(i,j,k ,PrimQ1_comp) : 0.0;
180  amrex::Real qc_hi = (n_qstate >= 2) ? cell_prim(i,j,k ,PrimQ2_comp) : 0.0;
181  amrex::Real qp_hi = (n_qstate >= 3) ? cell_prim(i,j,k ,PrimQ3_comp) : 0.0;
182 
183  amrex::Real qv_lo = (n_qstate >= 1) ? cell_prim(i,j,k-1,PrimQ1_comp) : 0.0;
184  amrex::Real qc_lo = (n_qstate >= 2) ? cell_prim(i,j,k-1,PrimQ2_comp) : 0.0;
185  amrex::Real qp_lo = (n_qstate >= 3) ? cell_prim(i,j,k-1,PrimQ3_comp) : 0.0;
186 
187  amrex::Real qv0_hi = qv0_arr(i,j,k );
188  amrex::Real qv0_lo = qv0_arr(i,j,k-1);
189 
190  amrex::Real q_hi = amrex::Real(0.61) * (qv_hi-qv0_hi) - (qc_hi + qp_hi)
191  + (cell_prim(i,j,k ,PrimTheta_comp) - th0_arr(i,j,k )) / th0_arr(i,j,k );
192  amrex::Real q_lo = amrex::Real(0.61) * (qv_lo-qv0_lo) - (qc_lo + qp_lo)
193  + (cell_prim(i,j,k-1,PrimTheta_comp) - th0_arr(i,j,k-1)) / th0_arr(i,j,k-1);
194  amrex::Real qavg = amrex::Real(0.5) * (q_hi + q_lo);
195 
196  amrex::Real r0avg = amrex::Real(0.5) * (r0_arr(i,j,k) + r0_arr(i,j,k-1));
197 
198  return ( -r0avg * grav_gpu * qavg);
199 }
200 
201 // **************************************************************************************
202 // Routines below this line are not currently used
203 // **************************************************************************************
204 
205 AMREX_GPU_DEVICE
206 AMREX_FORCE_INLINE
207 amrex::Real
208 buoyancy_dry_anelastic_T (int& i, int& j, int& k,
209  amrex::Real const& grav_gpu,
210  amrex::Real const& rd_over_cp,
211  const amrex::Array4<const amrex::Real>& r0_arr,
212  const amrex::Array4<const amrex::Real>& p0_arr,
213  const amrex::Array4<const amrex::Real>& cell_data)
214 {
215  amrex::Real rt0_hi = getRhoThetagivenP(p0_arr(i,j,k));
216  amrex::Real t0_hi = getTgivenPandTh(p0_arr(i,j,k), rt0_hi/r0_arr(i,j,k), rd_over_cp);
217  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);
218  amrex::Real q_hi = (t_hi-t0_hi)/t0_hi;
219 
220  amrex::Real rt0_lo = getRhoThetagivenP(p0_arr(i,j,k-1));
221  amrex::Real t0_lo = getTgivenPandTh(p0_arr(i,j,k-1), rt0_lo/r0_arr(i,j,k-1), rd_over_cp);
222  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);
223  amrex::Real q_lo = (t_lo-t0_lo)/t0_lo;
224 
225  amrex::Real r0_q_avg = amrex::Real(0.5) * (r0_arr(i,j,k) * q_hi + r0_arr(i,j,k-1) * q_lo);
226  return (-r0_q_avg * grav_gpu);
227 }
228 
229 #endif
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real buoyancy_rhopert(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 > &cell_data)
Definition: ERF_BuoyancyUtils.H:69
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:113
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)
Definition: ERF_BuoyancyUtils.H:171
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:87
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)
Definition: ERF_BuoyancyUtils.H:131
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)
Definition: ERF_BuoyancyUtils.H:30
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:208
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 PrimQ2_comp
Definition: ERF_IndexDefines.H:54
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define PrimQ3_comp
Definition: ERF_IndexDefines.H:55