Line data Source code
1 : !-------------------------------------------------------------------------
2 : !$Id$
3 : !===============================================================================
4 : module saturation
5 :
6 : ! Description:
7 : ! Contains functions that compute saturation with respect
8 : ! to liquid or ice.
9 : !-----------------------------------------------------------------------
10 :
11 : #ifdef GFDL
12 : use model_flags, only: & ! h1g, 2010-06-18
13 : I_sat_sphum
14 : #endif
15 :
16 : use clubb_precision, only: &
17 : core_rknd ! Variable(s)
18 :
19 :
20 : use model_flags, only: &
21 : saturation_formula, & ! Variable
22 : saturation_bolton, &
23 : saturation_gfdl, &
24 : saturation_flatau, &
25 : saturation_lookup
26 :
27 : implicit none
28 :
29 : private ! Change default so all items private
30 :
31 : public :: sat_mixrat_liq, sat_mixrat_ice, rcm_sat_adj, sat_vapor_press_liq
32 :
33 : private :: sat_vapor_press_ice_flatau, sat_vapor_press_ice_bolton
34 :
35 : interface sat_mixrat_liq
36 : module procedure sat_mixrat_liq_k ! Works over a single vertical level
37 : module procedure sat_mixrat_liq_2D ! Works over all vertical levels and columns
38 : end interface sat_mixrat_liq
39 :
40 : interface sat_mixrat_ice
41 : module procedure sat_mixrat_ice_k ! Works over a single vertical level
42 : module procedure sat_mixrat_ice_2D ! Works over all vertical levels and columns
43 : end interface sat_mixrat_ice
44 :
45 : ! Lookup table of values for saturation
46 : real( kind = core_rknd ), private, dimension(188:343) :: &
47 : svp_liq_lookup_table
48 :
49 : data svp_liq_lookup_table(188:343) / &
50 : 0.049560547_core_rknd, 0.059753418_core_rknd, 0.070129395_core_rknd, 0.083618164_core_rknd, &
51 : 0.09814453_core_rknd, 0.11444092_core_rknd, 0.13446045_core_rknd, 0.15686035_core_rknd, &
52 : 0.18218994_core_rknd, 0.21240234_core_rknd, 0.24725342_core_rknd, 0.28668213_core_rknd, &
53 : 0.33184814_core_rknd, 0.3826294_core_rknd, 0.4416504_core_rknd, 0.50775146_core_rknd, &
54 : 0.58343506_core_rknd, 0.6694946_core_rknd, 0.7668457_core_rknd, 0.87750244_core_rknd, &
55 : 1.0023804_core_rknd, 1.1434937_core_rknd, 1.3028564_core_rknd, 1.482544_core_rknd, &
56 : 1.6847534_core_rknd, 1.9118042_core_rknd, 2.1671143_core_rknd, 2.4535522_core_rknd, &
57 : 2.774231_core_rknd, 3.1330566_core_rknd, 3.5343628_core_rknd, 3.9819336_core_rknd, &
58 : 4.480713_core_rknd, 5.036072_core_rknd, 5.6540527_core_rknd, 6.340088_core_rknd, &
59 : 7.1015015_core_rknd, 7.9450684_core_rknd, 8.8793335_core_rknd, 9.91217_core_rknd, &
60 : 11.053528_core_rknd, 12.313049_core_rknd, 13.70166_core_rknd, 15.231018_core_rknd, &
61 : 16.91394_core_rknd, 18.764038_core_rknd, 20.795898_core_rknd, 23.025574_core_rknd, &
62 : 25.470093_core_rknd, 28.147766_core_rknd, 31.078003_core_rknd, 34.282043_core_rknd, &
63 : 37.782593_core_rknd, 41.60382_core_rknd, 45.771606_core_rknd, 50.31366_core_rknd, &
64 : 55.259644_core_rknd, 60.641174_core_rknd, 66.492004_core_rknd, 72.84802_core_rknd, &
65 : 79.74756_core_rknd, 87.23126_core_rknd, 95.34259_core_rknd, 104.12747_core_rknd, &
66 : 113.634796_core_rknd, 123.91641_core_rknd, 135.02725_core_rknd, 147.02563_core_rknd, &
67 : 159.97308_core_rknd, 173.93488_core_rknd, 188.97995_core_rknd, 205.18109_core_rknd, &
68 : 222.61517_core_rknd, 241.36334_core_rknd, 261.51108_core_rknd, 283.14853_core_rknd, &
69 : 306.37054_core_rknd, 331.27698_core_rknd, 357.97278_core_rknd, 386.56842_core_rknd, &
70 : 417.17978_core_rknd, 449.9286_core_rknd, 484.94254_core_rknd, 522.3556_core_rknd, &
71 : 562.30804_core_rknd, 604.947_core_rknd, 650.42645_core_rknd, 698.9074_core_rknd, &
72 : 750.55835_core_rknd, 805.55554_core_rknd, 864.0828_core_rknd, 926.3325_core_rknd, &
73 : 992.5052_core_rknd, 1062.8102_core_rknd, 1137.4657_core_rknd, 1216.6995_core_rknd, &
74 : 1300.7483_core_rknd, 1389.8594_core_rknd, 1484.2896_core_rknd, 1584.3064_core_rknd, &
75 : 1690.1881_core_rknd, 1802.224_core_rknd, 1920.7146_core_rknd, 2045.9724_core_rknd, &
76 : 2178.3218_core_rknd, 2318.099_core_rknd, 2465.654_core_rknd, 2621.3489_core_rknd, &
77 : 2785.5596_core_rknd, 2958.6758_core_rknd, 3141.101_core_rknd, 3333.2534_core_rknd, &
78 : 3535.5657_core_rknd, 3748.4863_core_rknd, 3972.4792_core_rknd, 4208.024_core_rknd, &
79 : 4455.616_core_rknd, 4715.7686_core_rknd, 4989.0127_core_rknd, 5275.8945_core_rknd, &
80 : 5576.9795_core_rknd, 5892.8535_core_rknd, 6224.116_core_rknd, 6571.3926_core_rknd, &
81 : 6935.3213_core_rknd, 7316.5674_core_rknd, 7715.8105_core_rknd, 8133.755_core_rknd, &
82 : 8571.125_core_rknd, 9028.667_core_rknd, 9507.15_core_rknd, 10007.367_core_rknd, &
83 : 10530.132_core_rknd, 11076.282_core_rknd, 11646.683_core_rknd, 12242.221_core_rknd, &
84 : 12863.808_core_rknd, 13512.384_core_rknd, 14188.913_core_rknd, 14894.385_core_rknd, &
85 : 15629.823_core_rknd, 16396.268_core_rknd, 17194.799_core_rknd, 18026.516_core_rknd, &
86 : 18892.55_core_rknd, 19794.07_core_rknd, 20732.262_core_rknd, 21708.352_core_rknd, &
87 : 22723.592_core_rknd, 23779.273_core_rknd, 24876.709_core_rknd, 26017.258_core_rknd, &
88 : 27202.3_core_rknd, 28433.256_core_rknd, 29711.578_core_rknd, 31038.766_core_rknd /
89 : !$acc declare create( svp_liq_lookup_table )
90 : !$omp threadprivate( svp_liq_lookup_table )
91 :
92 : contains
93 :
94 : !-------------------------------------------------------------------------
95 : ! Wrapped in interface sat_mixrat_liq
96 256360264 : function sat_mixrat_liq_k( p_in_Pa, T_in_K )
97 : !$acc routine seq
98 : ! Description:
99 : ! Used to compute the saturation mixing ratio of liquid water.
100 :
101 : ! References:
102 : ! Formula from Emanuel 1994, 4.4.14
103 : !-------------------------------------------------------------------------
104 :
105 : use constants_clubb, only: &
106 : ep
107 :
108 : use constants_clubb, only: T_freeze_K
109 :
110 : implicit none
111 :
112 : ! -------------------- Input Variables --------------------
113 : real( kind = core_rknd ), intent(in) :: &
114 : p_in_Pa, & ! Pressure [Pa]
115 : T_in_K ! Temperature [K]
116 :
117 : ! -------------------- Output Variables --------------------
118 : real( kind = core_rknd ) :: &
119 : sat_mixrat_liq_k
120 :
121 : ! -------------------- Local Variables --------------------
122 : real( kind = core_rknd ) :: &
123 : T_in_C, &
124 : T_in_C_sqd, &
125 : T_in_K_clipped
126 :
127 : integer :: &
128 : T_in_K_int
129 :
130 : ! Constant parameters
131 :
132 : ! Relative error norm expansion (-50 to 50 deg_C) from
133 : ! Table 3 of pp. 1510 of Flatau et al. 1992 (Water Vapor)
134 : ! (The 100 coefficient converts from mb to Pa)
135 : ! real, dimension(7), parameter :: a = &
136 : ! 100.* (/ 6.11176750, 0.443986062, 0.143053301E-01, &
137 : ! 0.265027242E-03, 0.302246994E-05, 0.203886313E-07, &
138 : ! 0.638780966E-10 /)
139 :
140 : ! Relative error norm expansion (-85 to 70 deg_C) from
141 : ! Table 4 of pp. 1511 of Flatau et al.
142 : !real( kind = core_rknd ), dimension(9), parameter :: a = &
143 : !100._core_rknd * &
144 : ! Commented out because the form has been redone, causing these number to no longer be needed,
145 : ! leaving them in for now for reference.
146 : ! (/ 6.11583699_core_rknd, 0.444606896_core_rknd, 0.143177157E-01_core_rknd, &
147 : ! 0.264224321E-03_core_rknd, 0.299291081E-05_core_rknd, 0.203154182E-07_core_rknd, &
148 : ! 0.702620698E-10_core_rknd, 0.379534310E-13_core_rknd,-0.321582393E-15_core_rknd /)
149 :
150 : real( kind = core_rknd ), parameter :: &
151 : min_T_in_C = -85._core_rknd, & ! [deg_C]
152 : min_T_in_K = 173.15_core_rknd ! Lowest temperature at which Goff-Gratch is valid [K]
153 :
154 : ! ---------------------- Output Variables ----------------------
155 : real( kind = core_rknd ) :: &
156 : esat ! Saturation vapor pressure over water [Pa]
157 :
158 : ! -------------------- Begin Code --------------------
159 :
160 : ! Calculate the SVP for water vapor.
161 0 : select case ( saturation_formula )
162 : case ( saturation_flatau )
163 :
164 : ! Using the Flatau, et al. polynomial approximation for SVP over vapor
165 :
166 : ! Determine deg K - 273.15
167 0 : T_in_C = T_in_K - T_freeze_K
168 :
169 : ! Since this approximation is only good out to -85 degrees Celsius we
170 : ! truncate the result here (Flatau, et al. 1992)
171 0 : T_in_C = max( T_in_C, min_T_in_C )
172 :
173 : ! Polynomial approx. (Flatau, et al. 1992)
174 :
175 : ! This is the generalized formula but is not computationally efficient.
176 : ! Based on Wexler's expressions(2.1)-(2.4) (See Flatau et al. p 1508)
177 : ! e_{sat} = a_1 + a_2 ( T - T_0 ) + ... + a_{n+1} ( T - T_0 )^n
178 :
179 : ! esat = a(1)
180 :
181 : ! do i = 2, size( a ) , 1
182 : ! esat = esat + a(i) * ( T_in_C )**(i-1)
183 : ! end do
184 :
185 : ! The 8th order polynomial fit. When running deep
186 : ! convective cases I noticed that absolute temperature often dips below
187 : ! -50 deg_C at higher altitudes, where the 6th order approximation is
188 : ! not accurate. -dschanen 20 Nov 2008
189 : !esat = a(1) + T_in_C*( a(2) + T_in_C*( a(3) + T_in_C*( a(4) + T_in_C &
190 : !*( a(5) + T_in_C*( a(6) + T_in_C*( a(7) + T_in_C*( a(8) + T_in_C*( a(9) ) ) ) ) ) ) ) )
191 :
192 :
193 : ! Factoring the polynomial above and changing it into this form allows the cpu
194 : ! to complete the calculations out of order. This is because modern cpus can complete
195 : ! multiple instructions at once if they do not depend on eachother, in the above case
196 : ! each instruction relies on the result of the last. In this version however, the terms
197 : ! in the parentheses could potentially be calculated in parallel by different execution
198 : ! units in the cpu, then only when those terms are being multiplied together do the
199 : ! instructions need to be done one at a time. See clubb issue 834 for more info.
200 : ! - Gunther Huebler, Aug 2018
201 0 : T_in_C_sqd = T_in_C**2
202 :
203 : esat = &
204 : - 3.21582393e-14_core_rknd * ( T_in_C - 646.5835252598777_core_rknd ) &
205 : * ( T_in_C + 90.72381630364440_core_rknd ) &
206 : * ( T_in_C_sqd + 111.0976961559954_core_rknd * T_in_C + 6459.629194243118_core_rknd ) &
207 : * ( T_in_C_sqd + 152.3131930092453_core_rknd * T_in_C + 6499.774954705265_core_rknd ) &
208 0 : * ( T_in_C_sqd + 174.4279584934021_core_rknd * T_in_C + 7721.679732114084_core_rknd )
209 :
210 : case ( saturation_bolton )
211 :
212 : ! Using the Bolton 1980 approximations for SVP over vapor
213 : ! Generally this more computationally expensive than the Flatau polnomial expansion
214 : esat = 611.2_core_rknd &
215 : * exp( (17.67_core_rknd *(T_in_K-T_freeze_K)) &
216 0 : / (T_in_K-29.65_core_rknd) ) ! Known magic number
217 :
218 : ! ---> h1g
219 : case ( saturation_gfdl )
220 :
221 : ! Using GFDL polynomial approximation for SVP with respect to liquid
222 : ! Since the Goff-Gratch approximation is valid only down to -70 degrees Celsius,
223 : ! we threshold the temperature. This will yield a minimal saturation at
224 : ! cold temperatures.
225 256360264 : T_in_K_clipped = max( min_T_in_K, T_in_K )
226 :
227 : ! Goff Gratch equation, uncertain below -70 C
228 :
229 : esat = 10._core_rknd**(-7.90298_core_rknd*(373.16_core_rknd/T_in_K_clipped-1._core_rknd)+ &
230 : 5.02808_core_rknd*log10(373.16_core_rknd/T_in_K_clipped)- &
231 : 1.3816e-7_core_rknd*(10._core_rknd**(11.344_core_rknd &
232 : *(1._core_rknd-T_in_K_clipped/373.16_core_rknd))-1._core_rknd)+ &
233 : 8.1328e-3_core_rknd*(10._core_rknd**(-3.49149_core_rknd &
234 : *(373.16_core_rknd/T_in_K_clipped-1._core_rknd))-1._core_rknd)+ &
235 256360264 : log10(1013.246_core_rknd))*100._core_rknd ! Known magic number
236 :
237 : ! <--- h1g
238 :
239 : case ( saturation_lookup )
240 :
241 0 : T_in_K_int = int( anint( T_in_K ) )
242 :
243 : ! Since this approximation is only good out to -85 degrees Celsius we
244 : ! truncate the result here
245 0 : T_in_K_int = min( max( T_in_K_int, 188 ), 343 )
246 :
247 : ! Use the lookup table to determine the saturation vapor pressure.
248 0 : esat = svp_liq_lookup_table( T_in_K_int )
249 :
250 : case default
251 :
252 : ! Undefined approximation
253 256360264 : esat = -99999.999_core_rknd
254 :
255 : end select
256 :
257 : ! If esat exceeds the air pressure, then assume esat~=0.5*pressure
258 : ! and set rsat = ep = 0.622
259 256360264 : if ( p_in_Pa-esat < 1.0_core_rknd ) then
260 : sat_mixrat_liq_k = ep
261 : else
262 :
263 : #ifdef GFDL
264 :
265 : ! GFDL uses specific humidity
266 : ! Formula for Saturation Specific Humidity
267 : if ( I_sat_sphum ) then ! h1g, 2010-06-18 begin mod
268 : sat_mixrat_liq_k = ep * ( esat / ( p_in_Pa &
269 : - (1.0_core_rknd-ep) * esat ) )
270 : else
271 : sat_mixrat_liq_k = ep * ( esat / ( p_in_Pa - esat ) )
272 : endif ! h1g, 2010-06-18 end mod
273 : #else
274 : ! Formula for Saturation Mixing Ratio:
275 : !
276 : ! rs = (epsilon) * [ esat / ( p - esat ) ];
277 : ! where epsilon = R_d / R_v
278 256360264 : sat_mixrat_liq_k = ep * esat / ( p_in_Pa - esat )
279 : #endif
280 : end if
281 :
282 : return
283 : end function sat_mixrat_liq_k
284 :
285 : !-------------------------------------------------------------------------
286 : !
287 3176496 : function sat_mixrat_liq_2D( nz, ngrdcol, p_in_Pa, T_in_K, &
288 : start_index_in )
289 : !
290 : ! Description:
291 : ! Used to compute the saturation mixing ratio of liquid water.
292 : !
293 : ! References:
294 : ! Formula from Emanuel 1994, 4.4.14
295 : !-------------------------------------------------------------------------
296 :
297 : use constants_clubb, only: &
298 : ep
299 :
300 : use constants_clubb, only: T_freeze_K
301 :
302 : implicit none
303 :
304 : ! -------------------- Input Variables --------------------
305 : integer, intent(in) :: &
306 : nz, &
307 : ngrdcol
308 :
309 : integer, intent(in), optional :: &
310 : start_index_in
311 :
312 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
313 : p_in_Pa, & ! Pressure [Pa]
314 : T_in_K ! Temperature [K]
315 :
316 : ! -------------------- Output Variables --------------------
317 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
318 : sat_mixrat_liq_2D
319 :
320 : ! -------------------- Local Variables --------------------
321 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
322 3176496 : esat
323 :
324 : real( kind = core_rknd ) :: &
325 : T_in_C, &
326 : T_in_C_sqd, &
327 : T_in_K_clipped
328 :
329 : integer :: &
330 : T_in_K_int
331 :
332 : integer :: i,k
333 :
334 : ! Constant parameters
335 :
336 : ! Relative error norm expansion (-50 to 50 deg_C) from
337 : ! Table 3 of pp. 1510 of Flatau et al. 1992 (Water Vapor)
338 : ! (The 100 coefficient converts from mb to Pa)
339 : ! real, dimension(7), parameter :: a = &
340 : ! 100.* (/ 6.11176750, 0.443986062, 0.143053301E-01, &
341 : ! 0.265027242E-03, 0.302246994E-05, 0.203886313E-07, &
342 : ! 0.638780966E-10 /)
343 :
344 : ! Relative error norm expansion (-85 to 70 deg_C) from
345 : ! Table 4 of pp. 1511 of Flatau et al.
346 : !real( kind = core_rknd ), dimension(9), parameter :: a = &
347 : !100._core_rknd * &
348 : ! Commented out because the form has been redone, causing these number to no longer be needed,
349 : ! leaving them in for now for reference.
350 : ! (/ 6.11583699_core_rknd, 0.444606896_core_rknd, 0.143177157E-01_core_rknd, &
351 : ! 0.264224321E-03_core_rknd, 0.299291081E-05_core_rknd, 0.203154182E-07_core_rknd, &
352 : ! 0.702620698E-10_core_rknd, 0.379534310E-13_core_rknd,-0.321582393E-15_core_rknd /)
353 :
354 : real( kind = core_rknd ), parameter :: &
355 : min_T_in_C = -85._core_rknd, & ! [deg_C]
356 : min_T_in_K = 173.15_core_rknd ! Lowest temperature at which Goff-Gratch is valid [K]
357 :
358 : integer :: &
359 : start_index
360 :
361 : ! -------------------- Begin Code --------------------
362 :
363 : !$acc data create(esat) &
364 : !$acc copyin(p_in_Pa,T_in_K), &
365 : !$acc copyout(sat_mixrat_liq_2D)
366 :
367 : ! start_index is an optional argument and
368 : ! used for choosing the sub-arrays
369 3176496 : if ( present(start_index_in) ) then
370 705888 : start_index = start_index_in
371 : else
372 : start_index = 1
373 : end if
374 :
375 : select case ( saturation_formula )
376 : case ( saturation_flatau )
377 :
378 : !$acc parallel loop gang vector collapse(2) default(present)
379 0 : do k = start_index, nz
380 0 : do i = 1, ngrdcol
381 :
382 : ! Determine deg K - 273.15
383 0 : T_in_C = T_in_K(i,k) - T_freeze_K
384 :
385 : ! Since this approximation is only good out to -85 degrees Celsius we
386 : ! truncate the result here (Flatau, et al. 1992)
387 0 : T_in_C = max( T_in_C, min_T_in_C )
388 :
389 : ! Polynomial approx. (Flatau, et al. 1992)
390 :
391 : ! This is the generalized formula but is not computationally efficient.
392 : ! Based on Wexler's expressions(2.1)-(2.4) (See Flatau et al. p 1508)
393 : ! e_{sat} = a_1 + a_2 ( T - T_0 ) + ... + a_{n+1} ( T - T_0 )^n
394 :
395 : ! esat = a(1)
396 :
397 : ! do i = 2, size( a ) , 1
398 : ! esat = esat + a(i) * ( T_in_C )**(i-1)
399 : ! end do
400 :
401 : ! The 8th order polynomial fit. When running deep
402 : ! convective cases I noticed that absolute temperature often dips below
403 : ! -50 deg_C at higher altitudes, where the 6th order approximation is
404 : ! not accurate. -dschanen 20 Nov 2008
405 : !esat = a(1) + T_in_C*( a(2) + T_in_C*( a(3) + T_in_C*( a(4) + T_in_C &
406 : !*( a(5) + T_in_C*( a(6) + T_in_C*( a(7) + T_in_C*( a(8) + T_in_C*( a(9) ) ) ) ) ) ) ) )
407 :
408 :
409 : ! Factoring the polynomial above and changing it into this form allows the cpu
410 : ! to complete the calculations out of order. This is because modern cpus can complete
411 : ! multiple instructions at once if they do not depend on eachother, in the above case
412 : ! each instruction relies on the result of the last. In this version however, the terms
413 : ! in the parentheses could potentially be calculated in parallel by different execution
414 : ! units in the cpu, then only when those terms are being multiplied together do the
415 : ! instructions need to be done one at a time. See clubb issue 834 for more info.
416 : ! - Gunther Huebler, Aug 2018
417 0 : T_in_C_sqd = T_in_C**2
418 :
419 : esat(i,k) = &
420 : - 3.21582393e-14_core_rknd * ( T_in_C - 646.5835252598777_core_rknd ) &
421 : * ( T_in_C + 90.72381630364440_core_rknd ) &
422 : * ( T_in_C_sqd + 111.0976961559954_core_rknd * T_in_C + 6459.629194243118_core_rknd ) &
423 : * ( T_in_C_sqd + 152.3131930092453_core_rknd * T_in_C + 6499.774954705265_core_rknd ) &
424 0 : * ( T_in_C_sqd + 174.4279584934021_core_rknd * T_in_C + 7721.679732114084_core_rknd )
425 : end do
426 : end do
427 : !$acc end parallel loop
428 :
429 : case ( saturation_bolton )
430 :
431 : ! Using the Bolton 1980 approximations for SVP over vapor
432 : ! Generally this more computationally expensive than the Flatau polnomial expansion
433 : !$acc parallel loop gang vector collapse(2) default(present)
434 0 : do k = start_index, nz
435 0 : do i = 1, ngrdcol
436 0 : esat(i,k) = 611.2_core_rknd &
437 : * exp( (17.67_core_rknd *(T_in_K(i,k)-T_freeze_K)) &
438 0 : / (T_in_K(i,k)-29.65_core_rknd) ) ! Known magic number
439 : end do
440 : end do
441 : !$acc end parallel loop
442 :
443 : ! ---> h1g
444 : case ( saturation_gfdl )
445 :
446 : ! Using GFDL polynomial approximation for SVP with respect to liquid
447 : !$acc parallel loop gang vector collapse(2) default(present)
448 272119824 : do k = start_index, nz
449 4493904624 : do i = 1, ngrdcol
450 :
451 : ! Since the Goff-Gratch approximation is valid only down to -70 degrees Celsius,
452 : ! we threshold the temperature. This will yield a minimal saturation at
453 : ! cold temperatures.
454 4221784800 : T_in_K_clipped = max( min_T_in_K, T_in_K(i,k) )
455 :
456 : ! Goff Gratch equation, uncertain below -70 C
457 :
458 : esat(i,k) = 10._core_rknd**(-7.90298_core_rknd*(373.16_core_rknd/T_in_K_clipped-1._core_rknd)+ &
459 : 5.02808_core_rknd*log10(373.16_core_rknd/T_in_K_clipped)- &
460 : 1.3816e-7_core_rknd*(10._core_rknd**(11.344_core_rknd &
461 : *(1._core_rknd-T_in_K_clipped/373.16_core_rknd))-1._core_rknd)+ &
462 : 8.1328e-3_core_rknd*(10._core_rknd**(-3.49149_core_rknd &
463 : *(373.16_core_rknd/T_in_K_clipped-1._core_rknd))-1._core_rknd)+ &
464 4490728128 : log10(1013.246_core_rknd))*100._core_rknd ! Known magic number
465 : end do
466 : end do
467 : !$acc end parallel loop
468 :
469 : ! <--- h1g
470 :
471 : case ( saturation_lookup )
472 :
473 : !$acc parallel loop gang vector collapse(2) default(present)
474 0 : do k = start_index, nz
475 0 : do i = 1, ngrdcol
476 0 : T_in_K_int = int( anint( T_in_K(i,k) ) )
477 :
478 : ! Since this approximation is only good out to -85 degrees Celsius we
479 : ! truncate the result here
480 0 : T_in_K_int = min( max( T_in_K_int, 188 ), 343 )
481 :
482 : ! Use the lookup table to determine the saturation vapor pressure.
483 0 : esat(i,k) = svp_liq_lookup_table( T_in_K_int )
484 : end do
485 : end do
486 : !$acc end parallel loop
487 :
488 : case default
489 :
490 : ! Undefined approximation
491 3176496 : esat = -99999.999_core_rknd
492 :
493 : end select
494 :
495 : !$acc parallel loop gang vector collapse(2) default(present)
496 272119824 : do k = start_index, nz
497 4493904624 : do i = 1, ngrdcol
498 :
499 : ! If esat exceeds the air pressure, then assume esat~=0.5*pressure
500 : ! and set rsat = ep = 0.622
501 4490728128 : if ( p_in_Pa(i,k)-esat(i,k) < 1.0_core_rknd ) then
502 80779423 : sat_mixrat_liq_2D(i,k) = ep
503 : else
504 :
505 : #ifdef GFDL
506 :
507 : ! GFDL uses specific humidity
508 : ! Formula for Saturation Specific Humidity
509 : if ( I_sat_sphum ) then ! h1g, 2010-06-18 begin mod
510 : sat_mixrat_liq_2D(i,k) = ep(i,k) * ( esat(i,k) / ( p_in_Pa(i,k) &
511 : - (1.0_core_rknd-ep) * esat(i,k) ) )
512 : else
513 : sat_mixrat_liq_2D(i,k) = ep(i,k) * ( esat(i,k) / ( p_in_Pa(i,k) - esat(i,k) ) )
514 : endif ! h1g, 2010-06-18 end mod
515 : #else
516 : ! Formula for Saturation Mixing Ratio:
517 : !
518 : ! rs = (epsilon) * [ esat / ( p - esat ) ];
519 : ! where epsilon = R_d / R_v
520 4141005377 : sat_mixrat_liq_2D(i,k) = ep * esat(i,k) / ( p_in_Pa(i,k) - esat(i,k) )
521 : #endif
522 : end if
523 :
524 : end do
525 : end do
526 : !$acc end parallel loop
527 :
528 : !$acc end data
529 :
530 : end function sat_mixrat_liq_2D
531 :
532 : !-----------------------------------------------------------------
533 : ! Wrapped in interface sat_vapor_press_liq
534 0 : subroutine sat_vapor_press_liq( T_in_K, &
535 : esat )
536 :
537 : ! Description:
538 : ! Computes SVP for water vapor. Calls one of the other functions
539 : ! that calculate an approximation to SVP.
540 :
541 : ! References:
542 : ! None
543 : !-----------------------------------------------------------------
544 : use model_flags, only: &
545 : saturation_formula, & ! Variable
546 : saturation_bolton, &
547 : saturation_gfdl, &
548 : saturation_flatau
549 :
550 : use clubb_precision, only: &
551 : core_rknd ! Variable(s)
552 :
553 : implicit none
554 :
555 : ! ------------------------ Input Variables ------------------------
556 : real( kind = core_rknd ), intent(in) :: &
557 : T_in_K ! Temperature [K]
558 :
559 : ! ------------------------ Output Variables ------------------------
560 : real( kind = core_rknd ), intent(out) :: &
561 : esat ! Saturation Vapor Pressure over Water [Pa]
562 :
563 : ! ------------------------ Being Code ------------------------
564 :
565 : ! Saturation Vapor Pressure, esat, can be found to be approximated
566 : ! in many different ways.
567 0 : select case ( saturation_formula )
568 : case ( saturation_bolton )
569 :
570 : ! Using the Bolton 1980 approximations for SVP over vapor
571 : call sat_vapor_press_liq_bolton( T_in_K, &
572 0 : esat )
573 :
574 : !Earthworks case
575 : case ( saturation_flatau )
576 :
577 : ! Using the Flatau, et al. polynomial approximation for SVP over vapor
578 : call sat_vapor_press_liq_flatau( T_in_K, &
579 0 : esat )
580 :
581 : ! ---> h1g
582 : case ( saturation_gfdl )
583 :
584 : ! Using GFDL polynomial approximation for SVP with respect to liquid
585 : call sat_vapor_press_liq_gfdl( T_in_K, &
586 0 : esat )
587 :
588 : ! <--- h1g
589 : case ( saturation_lookup )
590 :
591 : ! Use the lookup table to determine the saturation vapor pressure.
592 0 : esat = sat_vapor_press_liq_lookup( T_in_K )
593 :
594 : case default
595 :
596 : ! Undefined approximation
597 0 : esat = -99999.999_core_rknd
598 :
599 : end select
600 :
601 0 : return
602 :
603 : end subroutine sat_vapor_press_liq
604 :
605 : !------------------------------------------------------------------------
606 0 : subroutine sat_vapor_press_liq_flatau( T_in_K, &
607 : esat )
608 :
609 : ! Description:
610 : ! Computes SVP for water vapor.
611 :
612 : ! References:
613 : ! ``Polynomial Fits to Saturation Vapor Pressure'' Falatau, Walko,
614 : ! and Cotton. (1992) Journal of Applied Meteorology, Vol. 31,
615 : ! pp. 1507--1513
616 : !------------------------------------------------------------------------
617 :
618 : use constants_clubb, only: T_freeze_K
619 :
620 : use clubb_precision, only: &
621 : core_rknd ! Variable(s)
622 :
623 : implicit none
624 :
625 : ! Constant parameters
626 :
627 : ! Relative error norm expansion (-50 to 50 deg_C) from
628 : ! Table 3 of pp. 1510 of Flatau et al. 1992 (Water Vapor)
629 : ! (The 100 coefficient converts from mb to Pa)
630 : ! real, dimension(7), parameter :: a = &
631 : ! 100.* (/ 6.11176750, 0.443986062, 0.143053301E-01, &
632 : ! 0.265027242E-03, 0.302246994E-05, 0.203886313E-07, &
633 : ! 0.638780966E-10 /)
634 :
635 : ! Relative error norm expansion (-85 to 70 deg_C) from
636 : ! Table 4 of pp. 1511 of Flatau et al.
637 : !real( kind = core_rknd ), dimension(9), parameter :: a = &
638 : !100._core_rknd * &
639 : ! Commented out because the form has been redone, causing these number to no longer be needed,
640 : ! leaving them in for now for reference.
641 : ! (/ 6.11583699_core_rknd, 0.444606896_core_rknd, 0.143177157E-01_core_rknd, &
642 : ! 0.264224321E-03_core_rknd, 0.299291081E-05_core_rknd, 0.203154182E-07_core_rknd, &
643 : ! 0.702620698E-10_core_rknd, 0.379534310E-13_core_rknd,-0.321582393E-15_core_rknd /)
644 :
645 : real( kind = core_rknd ), parameter :: min_T_in_C = -85._core_rknd ! [deg_C]
646 :
647 : ! ---------------------- Input Variables ----------------------
648 : real( kind = core_rknd ), intent(in) :: &
649 : T_in_K ! Temperature [K]
650 :
651 : ! ---------------------- Output Variables ----------------------
652 : real( kind = core_rknd ), intent(out) :: &
653 : esat ! Saturation vapor pressure over water [Pa]
654 :
655 : ! ---------------------- Local Variables ----------------------
656 : real( kind = core_rknd ) :: T_in_C, T_in_C_sqd
657 :
658 : ! ---------------------- Begin Code ----------------------
659 :
660 : ! Determine deg K - 273.15
661 0 : T_in_C = T_in_K - T_freeze_K
662 :
663 : ! Since this approximation is only good out to -85 degrees Celsius we
664 : ! truncate the result here (Flatau, et al. 1992)
665 0 : T_in_C = max( T_in_C, min_T_in_C )
666 :
667 : ! Polynomial approx. (Flatau, et al. 1992)
668 :
669 : ! This is the generalized formula but is not computationally efficient.
670 : ! Based on Wexler's expressions(2.1)-(2.4) (See Flatau et al. p 1508)
671 : ! e_{sat} = a_1 + a_2 ( T - T_0 ) + ... + a_{n+1} ( T - T_0 )^n
672 :
673 : ! esat = a(1)
674 :
675 : ! do i = 2, size( a ) , 1
676 : ! esat = esat + a(i) * ( T_in_C )**(i-1)
677 : ! end do
678 :
679 : ! The 8th order polynomial fit. When running deep
680 : ! convective cases I noticed that absolute temperature often dips below
681 : ! -50 deg_C at higher altitudes, where the 6th order approximation is
682 : ! not accurate. -dschanen 20 Nov 2008
683 : !esat = a(1) + T_in_C*( a(2) + T_in_C*( a(3) + T_in_C*( a(4) + T_in_C &
684 : !*( a(5) + T_in_C*( a(6) + T_in_C*( a(7) + T_in_C*( a(8) + T_in_C*( a(9) ) ) ) ) ) ) ) )
685 :
686 :
687 : ! Factoring the polynomial above and changing it into this form allows the cpu
688 : ! to complete the calculations out of order. This is because modern cpus can complete
689 : ! multiple instructions at once if they do not depend on eachother, in the above case
690 : ! each instruction relies on the result of the last. In this version however, the terms
691 : ! in the parentheses could potentially be calculated in parallel by different execution
692 : ! units in the cpu, then only when those terms are being multiplied together do the
693 : ! instructions need to be done one at a time. See clubb issue 834 for more info.
694 : ! - Gunther Huebler, Aug 2018
695 0 : T_in_C_sqd = T_in_C**2
696 :
697 : esat = &
698 : - 3.21582393e-14_core_rknd * ( T_in_C - 646.5835252598777_core_rknd ) &
699 : * ( T_in_C + 90.72381630364440_core_rknd ) &
700 : * ( T_in_C_sqd + 111.0976961559954_core_rknd * T_in_C + 6459.629194243118_core_rknd ) &
701 : * ( T_in_C_sqd + 152.3131930092453_core_rknd * T_in_C + 6499.774954705265_core_rknd ) &
702 0 : * ( T_in_C_sqd + 174.4279584934021_core_rknd * T_in_C + 7721.679732114084_core_rknd )
703 :
704 0 : return
705 : end subroutine sat_vapor_press_liq_flatau
706 :
707 :
708 : !------------------------------------------------------------------------
709 0 : subroutine sat_vapor_press_liq_bolton( T_in_K, &
710 : esat )
711 : ! Description:
712 : ! Computes SVP for water vapor.
713 : ! References:
714 : ! Bolton 1980
715 : !------------------------------------------------------------------------
716 :
717 : use constants_clubb, only: T_freeze_K
718 :
719 : use clubb_precision, only: &
720 : core_rknd ! Variable(s)
721 :
722 : implicit none
723 :
724 : ! --------------------- Input Variables ---------------------
725 : real( kind = core_rknd ), intent(in) :: &
726 : T_in_K ! Temperature [K]
727 :
728 : ! --------------------- Output Variables ---------------------
729 : real( kind = core_rknd ), intent(out) :: &
730 : esat ! Saturation vapor pressure over water [Pa]
731 :
732 :
733 : ! --------------------------- Begin Code ---------------------------
734 :
735 : ! (Bolton 1980) approx.
736 : ! Generally this more computationally expensive than the Flatau polnomial expansion
737 : esat = 611.2_core_rknd &
738 : * exp( (17.67_core_rknd *(T_in_K-T_freeze_K)) &
739 0 : / (T_in_K-29.65_core_rknd) ) ! Known magic number
740 :
741 0 : return
742 : end subroutine sat_vapor_press_liq_bolton
743 :
744 :
745 : ! ---> h1g, 2010-06-16
746 : !------------------------------------------------------------------------
747 0 : subroutine sat_vapor_press_liq_gfdl( T_in_K, &
748 : esat )
749 : ! Description:
750 : ! copy from "GFDL polysvp.F90"
751 : ! Compute saturation vapor pressure with respect to liquid by using
752 : ! function from Goff and Gratch (1946)
753 :
754 : ! Polysvp returned in units of pa.
755 : ! T_in_K is input in units of K.
756 : !------------------------------------------------------------------------
757 :
758 : use clubb_precision, only: &
759 : core_rknd ! Variable(s)
760 :
761 : implicit none
762 :
763 : ! --------------------------- Input Variables ---------------------------
764 : real( kind = core_rknd ), intent(in) :: &
765 : T_in_K ! Absolute temperature [K]
766 :
767 : ! --------------------------- Output Variables ---------------------------
768 : real( kind = core_rknd ), intent(out) :: &
769 : esat ! Saturation vapor pressure over water [Pa]
770 :
771 : ! --------------------------- Local Variables ---------------------------
772 : real( kind = core_rknd ), parameter :: &
773 : min_T_in_K = 203.15_core_rknd ! Lowest temperature at which Goff-Gratch is valid [K]
774 :
775 : real( kind = core_rknd ) :: &
776 : T_in_K_clipped ! Absolute temperature with minimum threshold applied [K]
777 :
778 : ! --------------------------- Begin Code ---------------------------
779 :
780 : ! Since the Goff-Gratch approximation is valid only down to -70 degrees Celsius,
781 : ! we threshold the temperature. This will yield a minimal saturation at
782 : ! cold temperatures.
783 0 : T_in_K_clipped = max( min_T_in_K, T_in_K )
784 :
785 : ! Goff Gratch equation, uncertain below -70 C
786 :
787 : esat = 10._core_rknd**(-7.90298_core_rknd*(373.16_core_rknd/T_in_K_clipped-1._core_rknd)+ &
788 : 5.02808_core_rknd*log10(373.16_core_rknd/T_in_K_clipped)- &
789 : 1.3816e-7_core_rknd*(10._core_rknd**(11.344_core_rknd &
790 : *(1._core_rknd-T_in_K_clipped/373.16_core_rknd))-1._core_rknd)+ &
791 : 8.1328e-3_core_rknd*(10._core_rknd**(-3.49149_core_rknd &
792 : *(373.16_core_rknd/T_in_K_clipped-1._core_rknd))-1._core_rknd)+ &
793 0 : log10(1013.246_core_rknd))*100._core_rknd ! Known magic number
794 :
795 0 : return
796 : end subroutine sat_vapor_press_liq_gfdl
797 : ! <--- h1g, 2010-06-16
798 :
799 : !------------------------------------------------------------------------
800 0 : elemental function sat_vapor_press_liq_lookup( T_in_K ) result ( esat )
801 :
802 : ! Description:
803 : ! Computes SVP for water vapor, using a lookup table.
804 : !
805 : ! The lookup table was constructed using the Flatau approximation.
806 :
807 : ! References:
808 : ! ``Polynomial Fits to Saturation Vapor Pressure'' Falatau, Walko,
809 : ! and Cotton. (1992) Journal of Applied Meteorology, Vol. 31,
810 : ! pp. 1507--1513
811 : !------------------------------------------------------------------------
812 :
813 : implicit none
814 :
815 : ! External
816 : intrinsic :: max, min, int, anint
817 :
818 : ! Input Variables
819 : real( kind = core_rknd ), intent(in) :: T_in_K ! Temperature [K]
820 :
821 : ! Output Variables
822 : real( kind = core_rknd ) :: esat ! Saturation vapor pressure over water [Pa]
823 :
824 : ! Local Variables
825 : integer :: T_in_K_int
826 :
827 : ! ---- Begin Code ----
828 :
829 0 : T_in_K_int = int( anint( T_in_K ) )
830 :
831 : ! Since this approximation is only good out to -85 degrees Celsius we
832 : ! truncate the result here
833 0 : T_in_K_int = min( max( T_in_K_int, 188 ), 343 )
834 :
835 : ! Use the lookup table to determine the saturation vapor pressure.
836 0 : esat = svp_liq_lookup_table( T_in_K_int )
837 :
838 : return
839 : end function sat_vapor_press_liq_lookup
840 :
841 : !------------------------------------------------------------------------
842 : ! Wrapped in interface sat_mixrat_ice
843 0 : function sat_mixrat_ice_k( p_in_Pa, T_in_K )
844 :
845 : ! Description:
846 : ! Used to compute the saturation mixing ratio of ice.
847 :
848 : ! References:
849 : ! Formula from Emanuel 1994, 4.4.15
850 : !-------------------------------------------------------------------------
851 :
852 : use constants_clubb, only: &
853 : ep ! Variable(s)
854 :
855 : use clubb_precision, only: &
856 : core_rknd ! Variable(s)
857 :
858 : implicit none
859 :
860 : ! ------------------------ Input Variables ------------------------
861 : real( kind = core_rknd ), intent(in) :: &
862 : p_in_Pa, & ! Pressure [Pa]
863 : T_in_K ! Temperature [K]
864 :
865 : ! ------------------------ Output Variables ------------------------
866 : real( kind = core_rknd ) :: &
867 : sat_mixrat_ice_k
868 :
869 : ! ------------------------ Local Variables ------------------------
870 : real( kind = core_rknd ), dimension(1,1) :: &
871 : p_in_Pa_col, &
872 : T_in_K_col
873 :
874 : real( kind = core_rknd ), dimension(1,1) :: &
875 : sat_mixrat_ice_col
876 :
877 : integer :: i, k ! Loop indices
878 :
879 : ! ------------------------ Begin Code ------------------------
880 :
881 : ! Copy inputs to 2D arrays
882 0 : p_in_Pa_col(1,1) = p_in_Pa
883 0 : T_in_K_col(1,1) = T_in_K
884 :
885 : ! Call 2D version
886 0 : sat_mixrat_ice_col = sat_mixrat_ice_2D( 1, 1, p_in_Pa_col, T_in_K_col )
887 :
888 : ! Copy 2D result into output
889 0 : sat_mixrat_ice_k = sat_mixrat_ice_col(1,1)
890 :
891 : return
892 : end function sat_mixrat_ice_k
893 :
894 : !------------------------------------------------------------------------
895 : ! Wrapped in interface sat_mixrat_ice
896 1411776 : function sat_mixrat_ice_2D( nz, ngrdcol, p_in_Pa, T_in_K )
897 :
898 : ! Description:
899 : ! Used to compute the saturation mixing ratio of ice.
900 :
901 : ! References:
902 : ! Formula from Emanuel 1994, 4.4.15
903 : !-------------------------------------------------------------------------
904 :
905 : use constants_clubb, only: &
906 : T_freeze_K, & ! Variable(s)
907 : ep
908 :
909 : use clubb_precision, only: &
910 : core_rknd ! Variable(s)
911 :
912 : implicit none
913 :
914 : ! ------------------------ Input Variables ------------------------
915 : integer, intent(in) :: &
916 : nz, &
917 : ngrdcol
918 :
919 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
920 : p_in_Pa, & ! Pressure [Pa]
921 : T_in_K ! Temperature [K]
922 :
923 : ! ------------------------ Output Variables ------------------------
924 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
925 : sat_mixrat_ice_2D
926 :
927 : ! ------------------------ Local Variables ------------------------
928 : ! Relative error norm expansion (-90 to 0 deg_C) from
929 : ! Table 4 of pp. 1511 of Flatau et al. 1992 (Ice)
930 : real( kind = core_rknd ), dimension(9), parameter :: a = &
931 : 100._core_rknd * (/ 6.09868993_core_rknd, 0.499320233_core_rknd, 0.184672631E-01_core_rknd, &
932 : 0.402737184E-03_core_rknd, 0.565392987E-05_core_rknd, 0.521693933E-07_core_rknd, &
933 : 0.307839583E-09_core_rknd, 0.105785160E-11_core_rknd, 0.161444444E-14_core_rknd /)
934 :
935 : real( kind = core_rknd ), parameter :: &
936 : min_T_in_C = -90._core_rknd, & ! [deg_C]
937 : min_T_in_K = 173.15_core_rknd ! Lowest temperature at which Goff-Gratch is valid [K]
938 :
939 : real( kind = core_rknd ) :: &
940 : T_in_C, & ! Temperature [deg_C]
941 : T_in_K_clipped ! Absolute temperature with minimum threshold applied [K]
942 :
943 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
944 1411776 : esat_ice
945 :
946 : integer :: i, k ! Loop indices
947 :
948 : ! ------------------------ Begin Code ------------------------
949 :
950 : !$acc data create( esat_ice ) &
951 : !$acc copyin( p_in_Pa, T_in_K ) &
952 : !$acc copyout( sat_mixrat_ice_2D )
953 :
954 : ! Determine the SVP for the given temperature
955 : select case ( saturation_formula )
956 : case ( saturation_bolton )
957 :
958 : ! Using the Bolton 1980 approximations for SVP over ice
959 : !$acc parallel loop gang vector collapse(2) default(present)
960 0 : do i = 1, ngrdcol
961 0 : do k = 1, nz
962 :
963 : ! Exponential approx.
964 0 : esat_ice(i,k) = 100.0_core_rknd * exp( 23.33086_core_rknd - &
965 0 : (6111.72784_core_rknd/T_in_K(i,k)) + (0.15215_core_rknd*log( T_in_K(i,k) )) )
966 :
967 : end do
968 : end do
969 : !$acc end parallel loop
970 :
971 : case ( saturation_flatau )
972 :
973 : ! Using the Flatau, et al. polynomial approximation for SVP over ice
974 : !$acc parallel loop gang vector collapse(2) default(present)
975 0 : do i = 1, ngrdcol
976 0 : do k = 1, nz
977 :
978 : ! Determine deg K - 273.15
979 0 : T_in_C = T_in_K(i,k) - T_freeze_K
980 :
981 : ! Since this approximation is only good out to -90 degrees Celsius we
982 : ! truncate the result here (Flatau, et al. 1992)
983 0 : T_in_C = max( T_in_C, min_T_in_C )
984 :
985 : ! Polynomial approx. (Flatau, et al. 1992)
986 : ! esati = a(1)
987 :
988 : ! do i = 2, size( a ), 1
989 : ! esati = esati + a(i) * ( T_in_C )**(i-1)
990 : ! end do
991 :
992 : esat_ice(i,k) = a(1) + T_in_C*( a(2) + T_in_C*( a(3) + T_in_C*( a(4) + T_in_C &
993 0 : *( a(5) + T_in_C*( a(6) + T_in_C*( a(7) + T_in_C*( a(8) + T_in_C*a(9) ) ) ) ) ) ) )
994 :
995 : end do
996 : end do
997 : !$acc end parallel loop
998 :
999 : ! ---> h1g, 2010-06-16
1000 : case ( saturation_gfdl )
1001 :
1002 : ! Using GFDL polynomial approximation for SVP with respect to ice
1003 : !$acc parallel loop gang vector collapse(2) default(present)
1004 23573376 : do i = 1, ngrdcol
1005 1907309376 : do k = 1, nz
1006 :
1007 : ! Since the Goff-Gratch ice approximation is valid only down to -100 degrees Celsius,
1008 : ! we threshold the temperature. This will yield a minimal saturation at
1009 : ! cold temperatures.
1010 1883736000 : T_in_K_clipped = max( min_T_in_K, T_in_K(i,k) )
1011 :
1012 : ! Goff Gratch equation (good down to -100 C)
1013 :
1014 : esat_ice(i,k) = 10._core_rknd**(-9.09718_core_rknd* &
1015 : (273.16_core_rknd/T_in_K_clipped-1._core_rknd)-3.56654_core_rknd* &
1016 : log10(273.16_core_rknd/T_in_K_clipped)+0.876793_core_rknd* &
1017 : (1._core_rknd-T_in_K_clipped/273.16_core_rknd)+ &
1018 1905897600 : log10(6.1071_core_rknd))*100._core_rknd ! Known magic number
1019 : end do
1020 : end do
1021 : !$acc end parallel loop
1022 :
1023 : ! <--- h1g, 2010-06-16
1024 :
1025 : case default
1026 :
1027 : ! Undefined approximation
1028 1411776 : esat_ice = -99999.999_core_rknd
1029 :
1030 : end select
1031 :
1032 :
1033 : !$acc parallel loop gang vector collapse(2) default(present)
1034 121412736 : do k = 1, nz
1035 2005148736 : do i = 1, ngrdcol
1036 :
1037 : ! If esat_ice exceeds the air pressure, then assume esat_ice~=0.5*pressure
1038 : ! and set rsat = ep = 0.622
1039 2003736960 : if ( p_in_Pa(i,k)-esat_ice(i,k) < 1.0_core_rknd ) then
1040 35156389 : sat_mixrat_ice_2D(i,k) = ep
1041 : else
1042 :
1043 : #ifdef GFDL
1044 : ! GFDL uses specific humidity
1045 : ! Formula for Saturation Specific Humidity
1046 : if( I_sat_sphum ) then ! h1g, 2010-06-18 begin mod
1047 : sat_mixrat_ice_2D(i,k) = ep * ( esat_ice(i,k) &
1048 : / ( p_in_Pa(i,k) - (1.0_core_rknd-ep) * esat_ice(i,k) ) )
1049 : else
1050 : sat_mixrat_ice_2D(i,k) = ep * ( esat_ice(i,k) / ( p_in_Pa(i,k) - esat_ice(i,k) ) )
1051 : endif ! h1g, 2010-06-18 end mod
1052 : #else
1053 : ! Formula for Saturation Mixing Ratio:
1054 : !
1055 : ! rs = (epsilon) * [ esat / ( p - esat ) ];
1056 : ! where epsilon = R_d / R_v
1057 :
1058 1848579611 : sat_mixrat_ice_2D(i,k) = ep * esat_ice(i,k) / ( p_in_Pa(i,k) - esat_ice(i,k) )
1059 : #endif
1060 :
1061 : end if
1062 : end do
1063 : end do
1064 : !$acc end parallel loop
1065 :
1066 : !$acc end data
1067 :
1068 1411776 : return
1069 : end function sat_mixrat_ice_2D
1070 :
1071 : !------------------------------------------------------------------------
1072 : subroutine sat_vapor_press_ice( nz, ngrdcol, T_in_K, &
1073 : esat_ice )
1074 : !
1075 : ! Description:
1076 : ! Computes SVP for ice, using one of the various approximations.
1077 : !
1078 : ! References:
1079 : ! None
1080 : !------------------------------------------------------------------------
1081 :
1082 : use model_flags, only: &
1083 : saturation_formula, & ! Variable(s)
1084 : saturation_bolton, &
1085 : saturation_gfdl, &
1086 : saturation_flatau
1087 :
1088 : use clubb_precision, only: &
1089 : core_rknd ! Variable(s)
1090 :
1091 : implicit none
1092 :
1093 : ! ---------------------- Input Variable ----------------------
1094 : integer, intent(in) :: &
1095 : nz, &
1096 : ngrdcol
1097 :
1098 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
1099 : T_in_K ! Temperature [K]
1100 :
1101 : ! ---------------------- Output Variable ----------------------
1102 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
1103 : esat_ice ! Saturation Vapor Pressure over Ice [Pa]
1104 :
1105 : ! ---------------------- Begin Code ----------------------
1106 :
1107 : select case ( saturation_formula )
1108 : case ( saturation_bolton )
1109 :
1110 : ! Using the Bolton 1980 approximations for SVP over ice
1111 : call sat_vapor_press_ice_bolton( nz, ngrdcol, T_in_K, &
1112 : esat_ice )
1113 :
1114 : case ( saturation_flatau )
1115 :
1116 : ! Using the Flatau, et al. polynomial approximation for SVP over ice
1117 : call sat_vapor_press_ice_flatau( nz, ngrdcol, T_in_K, &
1118 : esat_ice )
1119 :
1120 : ! ---> h1g, 2010-06-16
1121 : case ( saturation_gfdl )
1122 :
1123 : ! Using GFDL polynomial approximation for SVP with respect to ice
1124 : call sat_vapor_press_ice_gfdl( nz, ngrdcol, T_in_K, &
1125 : esat_ice )
1126 : ! <--- h1g, 2010-06-16
1127 :
1128 : case default
1129 :
1130 : ! Undefined approximation
1131 : esat_ice = -99999.999_core_rknd
1132 :
1133 : end select
1134 :
1135 : return
1136 :
1137 : end subroutine sat_vapor_press_ice
1138 :
1139 : !------------------------------------------------------------------------
1140 : subroutine sat_vapor_press_ice_flatau( nz, ngrdcol, T_in_K, &
1141 : esati )
1142 : !
1143 : ! Description:
1144 : ! Computes SVP for ice.
1145 : !
1146 : ! References:
1147 : ! ``Polynomial Fits to Saturation Vapor Pressure'' Falatau, Walko,
1148 : ! and Cotton. (1992) Journal of Applied Meteorology, Vol. 31,
1149 : ! pp. 1507--1513
1150 : !------------------------------------------------------------------------
1151 : use constants_clubb, only: T_freeze_K
1152 :
1153 : use clubb_precision, only: &
1154 : core_rknd ! Variable(s)
1155 :
1156 : implicit none
1157 :
1158 : ! ------------------------ Input Variables ------------------------
1159 : integer, intent(in) :: &
1160 : nz, &
1161 : ngrdcol
1162 :
1163 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
1164 : T_in_K ! Temperature [deg_K]
1165 :
1166 : ! ------------------------ Output Variables ------------------------
1167 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
1168 : esati ! Saturation vapor pressure over ice [Pa]
1169 :
1170 : ! ------------------------ Local Variables ------------------------
1171 : ! Relative error norm expansion (-90 to 0 deg_C) from
1172 : ! Table 4 of pp. 1511 of Flatau et al. 1992 (Ice)
1173 : real( kind = core_rknd ), dimension(9), parameter :: a = &
1174 : 100._core_rknd * (/ 6.09868993_core_rknd, 0.499320233_core_rknd, 0.184672631E-01_core_rknd, &
1175 : 0.402737184E-03_core_rknd, 0.565392987E-05_core_rknd, 0.521693933E-07_core_rknd, &
1176 : 0.307839583E-09_core_rknd, 0.105785160E-11_core_rknd, 0.161444444E-14_core_rknd /)
1177 :
1178 : real( kind = core_rknd ), parameter :: min_T_in_C = -90._core_rknd ! [deg_C]
1179 :
1180 : real( kind = core_rknd ) :: &
1181 : T_in_C ! Temperature [deg_C]
1182 :
1183 : integer :: i, k
1184 :
1185 : ! ------------------------ Begin Code ------------------------
1186 :
1187 : do i = 1, ngrdcol
1188 : do k = 1, nz
1189 :
1190 : ! Determine deg K - 273.15
1191 : T_in_C = T_in_K(i,k) - T_freeze_K
1192 :
1193 : ! Since this approximation is only good out to -90 degrees Celsius we
1194 : ! truncate the result here (Flatau, et al. 1992)
1195 : T_in_C = max( T_in_C, min_T_in_C )
1196 :
1197 : ! Polynomial approx. (Flatau, et al. 1992)
1198 : ! esati = a(1)
1199 :
1200 : ! do i = 2, size( a ), 1
1201 : ! esati = esati + a(i) * ( T_in_C )**(i-1)
1202 : ! end do
1203 :
1204 : esati(i,k) = a(1) + T_in_C*( a(2) + T_in_C*( a(3) + T_in_C*( a(4) + T_in_C &
1205 : *( a(5) + T_in_C*( a(6) + T_in_C*( a(7) + T_in_C*( a(8) + T_in_C*a(9) ) ) ) ) ) ) )
1206 :
1207 : end do
1208 : end do
1209 :
1210 : return
1211 :
1212 : end subroutine sat_vapor_press_ice_flatau
1213 :
1214 : !------------------------------------------------------------------------
1215 : subroutine sat_vapor_press_ice_bolton( nz, ngrdcol, T_in_K, &
1216 : esati )
1217 : !
1218 : ! Description:
1219 : ! Computes SVP for ice.
1220 : !
1221 : ! References:
1222 : ! Bolton 1980
1223 : !------------------------------------------------------------------------
1224 :
1225 : use clubb_precision, only: &
1226 : core_rknd ! Variable(s)
1227 :
1228 : implicit none
1229 :
1230 : ! ------------------------ Input Variables ------------------------
1231 : integer, intent(in) :: &
1232 : nz, &
1233 : ngrdcol
1234 :
1235 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
1236 : T_in_K ! Temperature [deg_K]
1237 :
1238 : ! ------------------------ Output Variables ------------------------
1239 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
1240 : esati ! Saturation vapor pressure over ice [Pa]
1241 :
1242 : ! ------------------------ Local Variables ------------------------
1243 : integer :: i, k
1244 :
1245 : ! ------------------------ Begin Code ------------------------
1246 :
1247 : do i = 1, ngrdcol
1248 : do k = 1, nz
1249 :
1250 : ! Exponential approx.
1251 : esati(i,k) = 100.0_core_rknd * exp( 23.33086_core_rknd - &
1252 : (6111.72784_core_rknd/T_in_K(i,k)) + (0.15215_core_rknd*log( T_in_K(i,k) )) )
1253 :
1254 : end do
1255 : end do
1256 :
1257 : return
1258 :
1259 : end subroutine sat_vapor_press_ice_bolton
1260 :
1261 :
1262 : ! ---> h1g, 2010-06-16
1263 : !------------------------------------------------------------------------
1264 : subroutine sat_vapor_press_ice_gfdl( nz, ngrdcol, T_in_K, &
1265 : esati )
1266 : ! Description:
1267 : ! copy from "GFDL polysvp.F90"
1268 : ! Compute saturation vapor pressure with respect to liquid by using
1269 : ! function from Goff and Gratch (1946)
1270 : !
1271 : ! Polysvp returned in units of pa.
1272 : ! T_in_K is input in units of K.
1273 : !------------------------------------------------------------------------
1274 :
1275 : use clubb_precision, only: &
1276 : core_rknd ! Variable(s)
1277 :
1278 : implicit none
1279 :
1280 : ! ------------------------ Input Variables ------------------------
1281 : integer, intent(in) :: &
1282 : nz, &
1283 : ngrdcol
1284 :
1285 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
1286 : T_in_K ! Temperature [deg_K]
1287 :
1288 : ! ------------------------ Output Variables ------------------------
1289 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
1290 : esati ! Saturation vapor pressure over ice [Pa]
1291 :
1292 : ! ------------------------ Local Variables ------------------------
1293 : real( kind = core_rknd ), parameter :: &
1294 : min_T_in_K = 173.15_core_rknd ! Lowest temperature at which Goff-Gratch is valid [K]
1295 :
1296 : real( kind = core_rknd ) :: &
1297 : T_in_K_clipped ! Absolute temperature with minimum threshold applied [K]
1298 :
1299 : integer :: i, k
1300 :
1301 : ! ------------------------ Begin Code ------------------------
1302 :
1303 : do i = 1, ngrdcol
1304 : do k = 1, nz
1305 :
1306 : ! Since the Goff-Gratch ice approximation is valid only down to -100 degrees Celsius,
1307 : ! we threshold the temperature. This will yield a minimal saturation at
1308 : ! cold temperatures.
1309 : T_in_K_clipped = max( min_T_in_K, T_in_K(i,k) )
1310 :
1311 : ! Goff Gratch equation (good down to -100 C)
1312 :
1313 : esati(i,k) = 10._core_rknd**(-9.09718_core_rknd* &
1314 : (273.16_core_rknd/T_in_K_clipped-1._core_rknd)-3.56654_core_rknd* &
1315 : log10(273.16_core_rknd/T_in_K_clipped)+0.876793_core_rknd* &
1316 : (1._core_rknd-T_in_K_clipped/273.16_core_rknd)+ &
1317 : log10(6.1071_core_rknd))*100._core_rknd ! Known magic number
1318 : end do
1319 : end do
1320 :
1321 : return
1322 :
1323 : end subroutine sat_vapor_press_ice_gfdl
1324 : ! <--- h1g, 2010-06-16
1325 :
1326 : !-------------------------------------------------------------------------
1327 0 : function rcm_sat_adj( thlm, rtm, p_in_Pa, exner ) result ( rcm )
1328 :
1329 : ! Description:
1330 : !
1331 : ! This function uses an iterative method to find the value of rcm
1332 : ! from an initial profile that has saturation at some point.
1333 : !
1334 : ! References:
1335 : ! None
1336 : !-------------------------------------------------------------------------
1337 :
1338 : use clubb_precision, only: &
1339 : core_rknd ! Variable(s)
1340 :
1341 : use constants_clubb, only: &
1342 : Cp, & ! Variable(s)
1343 : Lv, &
1344 : zero_threshold
1345 :
1346 : implicit none
1347 :
1348 : ! Local Constant(s)
1349 : real( kind = core_rknd ), parameter :: &
1350 : tolerance = 0.001_core_rknd ! Tolerance on theta calculation [K]
1351 :
1352 : integer, parameter :: &
1353 : itermax = 1000000 ! Maximum interations
1354 :
1355 : ! External
1356 : intrinsic :: max, abs
1357 :
1358 : ! Input Variable(s)
1359 : real( kind = core_rknd ), intent(in) :: &
1360 : thlm, & ! Liquid Water Potential Temperature [K]
1361 : rtm, & ! Total Water Mixing Ratio [kg/kg]
1362 : p_in_Pa, & ! Pressure [Pa]
1363 : exner ! Exner function [-]
1364 :
1365 : ! Output Variable(s)
1366 : real( kind = core_rknd ) :: rcm ! Cloud water mixing ratio [kg/kg]
1367 :
1368 : ! Local Variable(s)
1369 : real( kind = core_rknd ) :: &
1370 : theta, answer, too_low, too_high, & ! [K]
1371 : rsat
1372 :
1373 : integer :: iteration
1374 :
1375 : ! ----- Begin Code -----
1376 :
1377 : ! Default initialization
1378 0 : theta = thlm
1379 0 : too_high = 0.0_core_rknd
1380 0 : too_low = 0.0_core_rknd
1381 :
1382 0 : do iteration = 1, itermax, 1
1383 :
1384 : answer = theta - (Lv/(Cp*exner)) &
1385 0 : *(MAX( rtm - sat_mixrat_liq(p_in_Pa,theta*exner), zero_threshold ))
1386 :
1387 0 : if ( ABS(answer - thlm) <= tolerance ) then
1388 : exit
1389 0 : else if ( answer - thlm > tolerance ) then
1390 : too_high = theta
1391 0 : else if ( thlm - answer > tolerance ) THEN
1392 0 : too_low = theta
1393 : end if
1394 :
1395 : ! For the first timestep, be sure to set a "too_high"
1396 : ! that is "way too high."
1397 0 : if ( iteration == 1 ) then
1398 0 : too_high = theta + 20.0_core_rknd
1399 : end if
1400 :
1401 0 : theta = (too_low + too_high)/2.0_core_rknd
1402 :
1403 : end do ! 1..itermax
1404 :
1405 0 : if ( iteration == itermax ) then
1406 : ! Magic Eric Raut added to remove compiler warning (clearly this value is not used)
1407 0 : rcm = 0.0_core_rknd
1408 :
1409 0 : error stop "Error in rcm_sat_adj: could not determine rcm"
1410 : else
1411 0 : rcm = MAX( rtm - sat_mixrat_liq( p_in_Pa, theta*exner), zero_threshold )
1412 : return
1413 : end if
1414 :
1415 : end function rcm_sat_adj
1416 :
1417 : end module saturation
|