-
Notifications
You must be signed in to change notification settings - Fork 0
/
wv_sat_methods.F90
485 lines (367 loc) · 13.6 KB
/
wv_sat_methods.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
module wv_sat_methods
! From: https://svn-ccsm-models.cgd.ucar.edu/cam1/trunk_tags/cam5_4_51/components/cam/src/physics/cam/wv_sat_methods.F90
! This portable module contains all CAM methods for estimating
! the saturation vapor pressure of water.
!
! wv_saturation provides CAM-specific interfaces and utilities
! based on these formulae.
!
! Typical usage of this module:
!
! Init:
! call wv_sat_methods_init(r8, <constants>, errstring)
!
! Get scheme index from a name string:
! scheme_idx = wv_sat_get_scheme_idx(scheme_name)
! if (.not. wv_sat_valid_idx(scheme_idx)) <throw some error>
!
! Get pressures:
! es = wv_sat_svp_water(t, scheme_idx)
! es = wv_sat_svp_ice(t, scheme_idx)
!
! Use ice/water transition range:
! es = wv_sat_svp_trice(t, ttrice, scheme_idx)
!
! Note that elemental functions cannot be pointed to, nor passed
! as arguments. If you need to do either, it is recommended to
! wrap the function so that it can be given an explicit (non-
! elemental) interface.
implicit none
private
save
integer, parameter :: r8 = selected_real_kind(12) ! 8 byte real
real(r8) :: tmelt ! Melting point of water at 1 atm (K)
real(r8) :: h2otrip ! Triple point temperature of water (K)
real(r8) :: tboil ! Boiling point of water at 1 atm (K)
real(r8) :: ttrice ! Ice-water transition range
real(r8) :: epsilo ! Ice-water transition range
real(r8) :: omeps ! 1._r8 - epsilo
! Indices representing individual schemes
integer, parameter :: Invalid_idx = -1
integer, parameter :: OldGoffGratch_idx = 0
integer, parameter :: GoffGratch_idx = 1
integer, parameter :: MurphyKoop_idx = 2
integer, parameter :: Bolton_idx = 3
! Index representing the current default scheme.
integer, parameter :: initial_default_idx = GoffGratch_idx
integer :: default_idx = initial_default_idx
public wv_sat_methods_init
public wv_sat_get_scheme_idx
public wv_sat_valid_idx
public wv_sat_set_default
public wv_sat_reset_default
public wv_sat_svp_water
public wv_sat_svp_ice
public wv_sat_svp_trans
! pressure -> humidity conversion
public wv_sat_svp_to_qsat
! Combined qsat operations
public wv_sat_qsat_water
public wv_sat_qsat_ice
public wv_sat_qsat_trans
contains
!---------------------------------------------------------------------
! ADMINISTRATIVE FUNCTIONS
!---------------------------------------------------------------------
! Get physical constants
subroutine wv_sat_methods_init(kind, tmelt_in, h2otrip_in, tboil_in, &
ttrice_in, epsilo_in, errstring)
integer, intent(in) :: kind
real(r8), intent(in) :: tmelt_in
real(r8), intent(in) :: h2otrip_in
real(r8), intent(in) :: tboil_in
real(r8), intent(in) :: ttrice_in
real(r8), intent(in) :: epsilo_in
character(len=*), intent(out) :: errstring
errstring = ' '
if (kind /= r8) then
write(errstring,*) 'wv_sat_methods_init: ERROR: ', &
kind,' was input kind but ',r8,' is internal kind.'
return
end if
if (ttrice_in < 0._r8) then
write(errstring,*) 'wv_sat_methods_init: ERROR: ', &
ttrice_in,' was input for ttrice, but negative range is invalid.'
return
end if
tmelt = tmelt_in
h2otrip = h2otrip_in
tboil = tboil_in
ttrice = ttrice_in
epsilo = epsilo_in
omeps = 1._r8 - epsilo
end subroutine wv_sat_methods_init
! Look up index by name.
pure function wv_sat_get_scheme_idx(name) result(idx)
character(len=*), intent(in) :: name
integer :: idx
select case (name)
case("GoffGratch")
idx = GoffGratch_idx
case("MurphyKoop")
idx = MurphyKoop_idx
case("OldGoffGratch")
idx = OldGoffGratch_idx
case("Bolton")
idx = Bolton_idx
case default
idx = Invalid_idx
end select
end function wv_sat_get_scheme_idx
! Check validity of an index from the above routine.
pure function wv_sat_valid_idx(idx) result(status)
integer, intent(in) :: idx
logical :: status
status = (idx /= Invalid_idx)
end function wv_sat_valid_idx
! Set default scheme (otherwise, Goff & Gratch is default)
! Returns a logical representing success (.true.) or
! failure (.false.).
function wv_sat_set_default(name) result(status)
character(len=*), intent(in) :: name
logical :: status
! Don't want to overwrite valid default with invalid,
! so assign to temporary and check it first.
integer :: tmp_idx
tmp_idx = wv_sat_get_scheme_idx(name)
status = wv_sat_valid_idx(tmp_idx)
if (status) default_idx = tmp_idx
end function wv_sat_set_default
! Reset default scheme to initial value.
! The same thing can be accomplished with wv_sat_set_default;
! the real reason to provide this routine is to reset the
! module for testing purposes.
subroutine wv_sat_reset_default()
default_idx = initial_default_idx
end subroutine wv_sat_reset_default
!---------------------------------------------------------------------
! UTILITIES
!---------------------------------------------------------------------
! Get saturation specific humidity given pressure and SVP.
! Specific humidity is limited to range 0-1.
elemental function wv_sat_svp_to_qsat(es, p) result(qs)
real(r8), intent(in) :: es ! SVP
real(r8), intent(in) :: p ! Current pressure.
real(r8) :: qs
! If pressure is less than SVP, set qs to maximum of 1.
if ( (p - es) <= 0._r8 ) then
qs = 1.0_r8
else
qs = epsilo*es / (p - omeps*es)
end if
end function wv_sat_svp_to_qsat
elemental subroutine wv_sat_qsat_water(t, p, es, qs, idx)
!------------------------------------------------------------------!
! Purpose: !
! Calculate SVP over water at a given temperature, and then !
! calculate and return saturation specific humidity. !
!------------------------------------------------------------------!
! Inputs
real(r8), intent(in) :: t ! Temperature
real(r8), intent(in) :: p ! Pressure
! Outputs
real(r8), intent(out) :: es ! Saturation vapor pressure
real(r8), intent(out) :: qs ! Saturation specific humidity
integer, intent(in), optional :: idx ! Scheme index
es = wv_sat_svp_water(t, idx)
qs = wv_sat_svp_to_qsat(es, p)
! Ensures returned es is consistent with limiters on qs.
es = min(es, p)
end subroutine wv_sat_qsat_water
elemental subroutine wv_sat_qsat_ice(t, p, es, qs, idx)
!------------------------------------------------------------------!
! Purpose: !
! Calculate SVP over ice at a given temperature, and then !
! calculate and return saturation specific humidity. !
!------------------------------------------------------------------!
! Inputs
real(r8), intent(in) :: t ! Temperature
real(r8), intent(in) :: p ! Pressure
! Outputs
real(r8), intent(out) :: es ! Saturation vapor pressure
real(r8), intent(out) :: qs ! Saturation specific humidity
integer, intent(in), optional :: idx ! Scheme index
es = wv_sat_svp_ice(t, idx)
qs = wv_sat_svp_to_qsat(es, p)
! Ensures returned es is consistent with limiters on qs.
es = min(es, p)
end subroutine wv_sat_qsat_ice
elemental subroutine wv_sat_qsat_trans(t, p, es, qs, idx)
!------------------------------------------------------------------!
! Purpose: !
! Calculate SVP over ice at a given temperature, and then !
! calculate and return saturation specific humidity. !
!------------------------------------------------------------------!
! Inputs
real(r8), intent(in) :: t ! Temperature
real(r8), intent(in) :: p ! Pressure
! Outputs
real(r8), intent(out) :: es ! Saturation vapor pressure
real(r8), intent(out) :: qs ! Saturation specific humidity
integer, intent(in), optional :: idx ! Scheme index
es = wv_sat_svp_trans(t, idx)
qs = wv_sat_svp_to_qsat(es, p)
! Ensures returned es is consistent with limiters on qs.
es = min(es, p)
end subroutine wv_sat_qsat_trans
!---------------------------------------------------------------------
! SVP INTERFACE FUNCTIONS
!---------------------------------------------------------------------
elemental function wv_sat_svp_water(t, idx) result(es)
real(r8), intent(in) :: t
integer, intent(in), optional :: idx
real(r8) :: es
integer :: use_idx
if (present(idx)) then
use_idx = idx
else
use_idx = default_idx
end if
select case (use_idx)
case(GoffGratch_idx)
es = GoffGratch_svp_water(t)
case(MurphyKoop_idx)
es = MurphyKoop_svp_water(t)
case(OldGoffGratch_idx)
es = OldGoffGratch_svp_water(t)
case(Bolton_idx)
es = Bolton_svp_water(t)
end select
end function wv_sat_svp_water
elemental function wv_sat_svp_ice(t, idx) result(es)
real(r8), intent(in) :: t
integer, intent(in), optional :: idx
real(r8) :: es
integer :: use_idx
if (present(idx)) then
use_idx = idx
else
use_idx = default_idx
end if
select case (use_idx)
case(GoffGratch_idx)
es = GoffGratch_svp_ice(t)
case(MurphyKoop_idx)
es = MurphyKoop_svp_ice(t)
case(OldGoffGratch_idx)
es = OldGoffGratch_svp_ice(t)
case(Bolton_idx)
es = Bolton_svp_water(t)
end select
end function wv_sat_svp_ice
elemental function wv_sat_svp_trans(t, idx) result (es)
real(r8), intent(in) :: t
integer, intent(in), optional :: idx
real(r8) :: es
real(r8) :: esice ! Saturation vapor pressure over ice
real(r8) :: weight ! Intermediate scratch variable for es transition
!
! Water
!
if (t >= (tmelt - ttrice)) then
es = wv_sat_svp_water(t,idx)
else
es = 0.0_r8
end if
!
! Ice
!
if (t < tmelt) then
esice = wv_sat_svp_ice(t,idx)
if ( (tmelt - t) > ttrice ) then
weight = 1.0_r8
else
weight = (tmelt - t)/ttrice
end if
es = weight*esice + (1.0_r8 - weight)*es
end if
end function wv_sat_svp_trans
!---------------------------------------------------------------------
! SVP METHODS
!---------------------------------------------------------------------
! Goff & Gratch (1946)
elemental function GoffGratch_svp_water(t) result(es)
real(r8), intent(in) :: t ! Temperature in Kelvin
real(r8) :: es ! SVP in Pa
! uncertain below -70 C
es = 10._r8**(-7.90298_r8*(tboil/t-1._r8)+ &
5.02808_r8*log10(tboil/t)- &
1.3816e-7_r8*(10._r8**(11.344_r8*(1._r8-t/tboil))-1._r8)+ &
8.1328e-3_r8*(10._r8**(-3.49149_r8*(tboil/t-1._r8))-1._r8)+ &
log10(1013.246_r8))*100._r8
end function GoffGratch_svp_water
elemental function GoffGratch_svp_ice(t) result(es)
real(r8), intent(in) :: t ! Temperature in Kelvin
real(r8) :: es ! SVP in Pa
! good down to -100 C
es = 10._r8**(-9.09718_r8*(h2otrip/t-1._r8)-3.56654_r8* &
log10(h2otrip/t)+0.876793_r8*(1._r8-t/h2otrip)+ &
log10(6.1071_r8))*100._r8
end function GoffGratch_svp_ice
! Murphy & Koop (2005)
elemental function MurphyKoop_svp_water(t) result(es)
real(r8), intent(in) :: t ! Temperature in Kelvin
real(r8) :: es ! SVP in Pa
! (good for 123 < T < 332 K)
es = exp(54.842763_r8 - (6763.22_r8 / t) - (4.210_r8 * log(t)) + &
(0.000367_r8 * t) + (tanh(0.0415_r8 * (t - 218.8_r8)) * &
(53.878_r8 - (1331.22_r8 / t) - (9.44523_r8 * log(t)) + &
0.014025_r8 * t)))
end function MurphyKoop_svp_water
elemental function MurphyKoop_svp_ice(t) result(es)
real(r8), intent(in) :: t ! Temperature in Kelvin
real(r8) :: es ! SVP in Pa
! (good down to 110 K)
es = exp(9.550426_r8 - (5723.265_r8 / t) + (3.53068_r8 * log(t)) &
- (0.00728332_r8 * t))
end function MurphyKoop_svp_ice
! Old CAM implementation, also labelled Goff & Gratch (1946)
! The water formula differs only due to compiler-dependent order of
! operations, so differences are roundoff level, usually 0.
! The ice formula gives fairly close answers to the current
! implementation, but has been rearranged, and uses the
! 1 atm melting point of water as the triple point.
! Differences are thus small but above roundoff.
! A curious fact: although using the melting point of water was
! probably a mistake, it mildly improves accuracy for ice svp,
! since it compensates for a systematic error in Goff & Gratch.
elemental function OldGoffGratch_svp_water(t) result(es)
real(r8), intent(in) :: t
real(r8) :: es
real(r8) :: ps, e1, e2, f1, f2, f3, f4, f5, f
ps = 1013.246_r8
e1 = 11.344_r8*(1.0_r8 - t/tboil)
e2 = -3.49149_r8*(tboil/t - 1.0_r8)
f1 = -7.90298_r8*(tboil/t - 1.0_r8)
f2 = 5.02808_r8*log10(tboil/t)
f3 = -1.3816_r8*(10.0_r8**e1 - 1.0_r8)/10000000.0_r8
f4 = 8.1328_r8*(10.0_r8**e2 - 1.0_r8)/1000.0_r8
f5 = log10(ps)
f = f1 + f2 + f3 + f4 + f5
es = (10.0_r8**f)*100.0_r8
end function OldGoffGratch_svp_water
elemental function OldGoffGratch_svp_ice(t) result(es)
real(r8), intent(in) :: t
real(r8) :: es
real(r8) :: term1, term2, term3
term1 = 2.01889049_r8/(tmelt/t)
term2 = 3.56654_r8*log(tmelt/t)
term3 = 20.947031_r8*(tmelt/t)
es = 575.185606e10_r8*exp(-(term1 + term2 + term3))
end function OldGoffGratch_svp_ice
! Bolton (1980)
! zm_conv deep convection scheme contained this SVP calculation.
! It appears to be from D. Bolton, 1980, Monthly Weather Review.
! Unlike the other schemes, no distinct ice formula is associated
! with it. (However, a Bolton ice formula exists in CLUBB.)
! The original formula used degrees C, but this function
! takes Kelvin and internally converts.
elemental function Bolton_svp_water(t) result(es)
real(r8),parameter :: c1 = 611.2_r8
real(r8),parameter :: c2 = 17.67_r8
real(r8),parameter :: c3 = 243.5_r8
real(r8), intent(in) :: t ! Temperature in Kelvin
real(r8) :: es ! SVP in Pa
es = c1*exp( (c2*(t - tmelt))/((t - tmelt)+c3) )
end function Bolton_svp_water
end module wv_sat_methods