1 | #define WRF_PORT |
---|
2 | |
---|
3 | ! |
---|
4 | ! Common block and statement functions for saturation vapor pressure |
---|
5 | ! look-up procedure, J. J. Hack, February 1990 |
---|
6 | ! |
---|
7 | ! Ported to WRF by William.Gustafson@pnl.gov, Nov. 2009 |
---|
8 | ! Updated to version from CESM1_0_1, Nov. 2010 |
---|
9 | ! |
---|
10 | ! $Id$ |
---|
11 | ! |
---|
12 | module wv_saturation |
---|
13 | use shr_kind_mod, only: r8 => shr_kind_r8 |
---|
14 | #ifndef WRF_PORT |
---|
15 | use abortutils, only: endrun |
---|
16 | use cam_logfile, only: iulog |
---|
17 | #else |
---|
18 | use module_cam_support, only: & |
---|
19 | masterproc, & |
---|
20 | endrun, & |
---|
21 | iulog |
---|
22 | use module_wrf_error |
---|
23 | #endif |
---|
24 | |
---|
25 | implicit none |
---|
26 | private |
---|
27 | save |
---|
28 | ! |
---|
29 | ! Public interfaces |
---|
30 | ! |
---|
31 | public gestbl ! Initialization subroutine |
---|
32 | public estblf ! saturation pressure table lookup |
---|
33 | public aqsat ! Returns saturation vapor pressure |
---|
34 | #ifndef WRF_PORT |
---|
35 | public aqsatd ! Same as aqsat, but also returns a temperature derivitive |
---|
36 | public vqsatd ! Vector version of aqsatd |
---|
37 | #endif |
---|
38 | public fqsatd ! Function version of vqsatd |
---|
39 | #ifndef WRF_PORT |
---|
40 | public qsat_water ! saturation mixing ration with respect to liquid water |
---|
41 | public vqsat_water ! vector version of qsat_water |
---|
42 | public qsat_ice ! saturation mixing ration with respect to ice |
---|
43 | public vqsat_ice ! vector version of qsat_ice |
---|
44 | public vqsatd_water |
---|
45 | public aqsat_water |
---|
46 | public vqsatd2_water ! Variant of vqsatd_water to print out dqsdT |
---|
47 | public vqsatd2_water_single ! Single value version of vqsatd2_water |
---|
48 | public vqsatd2 |
---|
49 | public vqsatd2_single |
---|
50 | public polysvp |
---|
51 | #endif |
---|
52 | ! |
---|
53 | ! Data used by cldwat |
---|
54 | ! |
---|
55 | public hlatv, tmin, hlatf, rgasv, pcf, cp, epsqs, ttrice |
---|
56 | ! |
---|
57 | ! Data |
---|
58 | ! |
---|
59 | integer plenest ! length of saturation vapor pressure table |
---|
60 | parameter (plenest=250) |
---|
61 | ! |
---|
62 | ! Table of saturation vapor pressure values es from tmin degrees |
---|
63 | ! to tmax+1 degrees k in one degree increments. ttrice defines the |
---|
64 | ! transition region where es is a combination of ice & water values |
---|
65 | ! |
---|
66 | real(r8) estbl(plenest) ! table values of saturation vapor pressure |
---|
67 | real(r8) tmin ! min temperature (K) for table |
---|
68 | real(r8) tmax ! max temperature (K) for table |
---|
69 | real(r8) ttrice ! transition range from es over H2O to es over ice |
---|
70 | real(r8) pcf(6) ! polynomial coeffs -> es transition water to ice |
---|
71 | real(r8) epsqs ! Ratio of h2o to dry air molecular weights |
---|
72 | real(r8) rgasv ! Gas constant for water vapor |
---|
73 | real(r8) hlatf ! Latent heat of vaporization |
---|
74 | real(r8) hlatv ! Latent heat of fusion |
---|
75 | real(r8) cp ! specific heat of dry air |
---|
76 | real(r8) tmelt ! Melting point of water (K) |
---|
77 | logical icephs ! false => saturation vapor press over water only |
---|
78 | |
---|
79 | contains |
---|
80 | |
---|
81 | real(r8) function estblf( td ) |
---|
82 | ! |
---|
83 | ! Saturation vapor pressure table lookup |
---|
84 | ! |
---|
85 | real(r8), intent(in) :: td ! Temperature for saturation lookup |
---|
86 | ! |
---|
87 | real(r8) :: e ! intermediate variable for es look-up |
---|
88 | real(r8) :: ai |
---|
89 | integer :: i |
---|
90 | ! |
---|
91 | e = max(min(td,tmax),tmin) ! partial pressure |
---|
92 | i = int(e-tmin)+1 |
---|
93 | ai = aint(e-tmin) |
---|
94 | estblf = (tmin+ai-e+1._r8)* & |
---|
95 | estbl(i)-(tmin+ai-e)* & |
---|
96 | estbl(i+1) |
---|
97 | end function estblf |
---|
98 | |
---|
99 | subroutine gestbl(tmn ,tmx ,trice ,ip ,epsil , & |
---|
100 | latvap ,latice ,rh2o ,cpair ,tmeltx ) |
---|
101 | !----------------------------------------------------------------------- |
---|
102 | ! |
---|
103 | ! Purpose: |
---|
104 | ! Builds saturation vapor pressure table for later lookup procedure. |
---|
105 | ! |
---|
106 | ! Method: |
---|
107 | ! Uses Goff & Gratch (1946) relationships to generate the table |
---|
108 | ! according to a set of free parameters defined below. Auxiliary |
---|
109 | ! routines are also included for making rapid estimates (well with 1%) |
---|
110 | ! of both es and d(es)/dt for the particular table configuration. |
---|
111 | ! |
---|
112 | ! Author: J. Hack |
---|
113 | ! |
---|
114 | !----------------------------------------------------------------------- |
---|
115 | #ifndef WRF_PORT |
---|
116 | use spmd_utils, only: masterproc |
---|
117 | #else |
---|
118 | use module_cam_gffgch |
---|
119 | #endif |
---|
120 | !------------------------------Arguments-------------------------------- |
---|
121 | ! |
---|
122 | ! Input arguments |
---|
123 | ! |
---|
124 | real(r8), intent(in) :: tmn ! Minimum temperature entry in es lookup table |
---|
125 | real(r8), intent(in) :: tmx ! Maximum temperature entry in es lookup table |
---|
126 | real(r8), intent(in) :: epsil ! Ratio of h2o to dry air molecular weights |
---|
127 | real(r8), intent(in) :: trice ! Transition range from es over range to es over ice |
---|
128 | real(r8), intent(in) :: latvap ! Latent heat of vaporization |
---|
129 | real(r8), intent(in) :: latice ! Latent heat of fusion |
---|
130 | real(r8), intent(in) :: rh2o ! Gas constant for water vapor |
---|
131 | real(r8), intent(in) :: cpair ! Specific heat of dry air |
---|
132 | real(r8), intent(in) :: tmeltx ! Melting point of water (K) |
---|
133 | ! |
---|
134 | !---------------------------Local variables----------------------------- |
---|
135 | ! |
---|
136 | real(r8) t ! Temperature |
---|
137 | integer n ! Increment counter |
---|
138 | integer lentbl ! Calculated length of lookup table |
---|
139 | integer itype ! Ice phase: 0 -> no ice phase |
---|
140 | ! 1 -> ice phase, no transition |
---|
141 | ! -x -> ice phase, x degree transition |
---|
142 | logical ip ! Ice phase logical flag |
---|
143 | ! |
---|
144 | !----------------------------------------------------------------------- |
---|
145 | ! |
---|
146 | ! Set es table parameters |
---|
147 | ! |
---|
148 | tmin = tmn ! Minimum temperature entry in table |
---|
149 | tmax = tmx ! Maximum temperature entry in table |
---|
150 | ttrice = trice ! Trans. range from es over h2o to es over ice |
---|
151 | icephs = ip ! Ice phase (true or false) |
---|
152 | ! |
---|
153 | ! Set physical constants required for es calculation |
---|
154 | ! |
---|
155 | epsqs = epsil |
---|
156 | hlatv = latvap |
---|
157 | hlatf = latice |
---|
158 | rgasv = rh2o |
---|
159 | cp = cpair |
---|
160 | tmelt = tmeltx |
---|
161 | ! |
---|
162 | lentbl = INT(tmax-tmin+2.000001_r8) |
---|
163 | if (lentbl .gt. plenest) then |
---|
164 | write(iulog,9000) tmax, tmin, plenest |
---|
165 | #ifdef WRF_PORT |
---|
166 | call wrf_message(iulog) |
---|
167 | #endif |
---|
168 | call endrun ('GESTBL') ! Abnormal termination |
---|
169 | end if |
---|
170 | ! |
---|
171 | ! Begin building es table. |
---|
172 | ! Check whether ice phase requested. |
---|
173 | ! If so, set appropriate transition range for temperature |
---|
174 | ! |
---|
175 | if (icephs) then |
---|
176 | if (ttrice /= 0.0_r8) then |
---|
177 | itype = -ttrice |
---|
178 | else |
---|
179 | itype = 1 |
---|
180 | end if |
---|
181 | else |
---|
182 | itype = 0 |
---|
183 | end if |
---|
184 | ! |
---|
185 | t = tmin - 1.0_r8 |
---|
186 | do n=1,lentbl |
---|
187 | t = t + 1.0_r8 |
---|
188 | call gffgch(t,estbl(n),itype) |
---|
189 | end do |
---|
190 | ! |
---|
191 | do n=lentbl+1,plenest |
---|
192 | estbl(n) = -99999.0_r8 |
---|
193 | end do |
---|
194 | ! |
---|
195 | ! Table complete -- Set coefficients for polynomial approximation of |
---|
196 | ! difference between saturation vapor press over water and saturation |
---|
197 | ! pressure over ice for -ttrice < t < 0 (degrees C). NOTE: polynomial |
---|
198 | ! is valid in the range -40 < t < 0 (degrees C). |
---|
199 | ! |
---|
200 | ! --- Degree 5 approximation --- |
---|
201 | ! |
---|
202 | pcf(1) = 5.04469588506e-01_r8 |
---|
203 | pcf(2) = -5.47288442819e+00_r8 |
---|
204 | pcf(3) = -3.67471858735e-01_r8 |
---|
205 | pcf(4) = -8.95963532403e-03_r8 |
---|
206 | pcf(5) = -7.78053686625e-05_r8 |
---|
207 | ! |
---|
208 | ! --- Degree 6 approximation --- |
---|
209 | ! |
---|
210 | !-----pcf(1) = 7.63285250063e-02 |
---|
211 | !-----pcf(2) = -5.86048427932e+00 |
---|
212 | !-----pcf(3) = -4.38660831780e-01 |
---|
213 | !-----pcf(4) = -1.37898276415e-02 |
---|
214 | !-----pcf(5) = -2.14444472424e-04 |
---|
215 | !-----pcf(6) = -1.36639103771e-06 |
---|
216 | ! |
---|
217 | if (masterproc) then |
---|
218 | write(iulog,*)' *** SATURATION VAPOR PRESSURE TABLE COMPLETED ***' |
---|
219 | #ifdef WRF_PORT |
---|
220 | call wrf_debug(100,iulog) |
---|
221 | #endif |
---|
222 | end if |
---|
223 | |
---|
224 | return |
---|
225 | ! |
---|
226 | 9000 format('GESTBL: FATAL ERROR *********************************',/, & |
---|
227 | ' TMAX AND TMIN REQUIRE A LARGER DIMENSION ON THE LENGTH', & |
---|
228 | ' OF THE SATURATION VAPOR PRESSURE TABLE ESTBL(PLENEST)',/, & |
---|
229 | ' TMAX, TMIN, AND PLENEST => ', 2f7.2, i3) |
---|
230 | ! |
---|
231 | end subroutine gestbl |
---|
232 | |
---|
233 | subroutine aqsat(t ,p ,es ,qs ,ii , & |
---|
234 | ilen ,kk ,kstart ,kend ) |
---|
235 | !----------------------------------------------------------------------- |
---|
236 | ! |
---|
237 | ! Purpose: |
---|
238 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
239 | ! precomputed table, calculate and return saturation specific humidity |
---|
240 | ! (g/g),for input arrays of temperature and pressure (dimensioned ii,kk) |
---|
241 | ! This routine is useful for evaluating only a selected region in the |
---|
242 | ! vertical. |
---|
243 | ! |
---|
244 | ! Method: |
---|
245 | ! <Describe the algorithm(s) used in the routine.> |
---|
246 | ! <Also include any applicable external references.> |
---|
247 | ! |
---|
248 | ! Author: J. Hack |
---|
249 | ! |
---|
250 | !------------------------------Arguments-------------------------------- |
---|
251 | ! |
---|
252 | ! Input arguments |
---|
253 | ! |
---|
254 | integer, intent(in) :: ii ! I dimension of arrays t, p, es, qs |
---|
255 | integer, intent(in) :: kk ! K dimension of arrays t, p, es, qs |
---|
256 | integer, intent(in) :: ilen ! Length of vectors in I direction which |
---|
257 | integer, intent(in) :: kstart ! Starting location in K direction |
---|
258 | integer, intent(in) :: kend ! Ending location in K direction |
---|
259 | real(r8), intent(in) :: t(ii,kk) ! Temperature |
---|
260 | real(r8), intent(in) :: p(ii,kk) ! Pressure |
---|
261 | ! |
---|
262 | ! Output arguments |
---|
263 | ! |
---|
264 | real(r8), intent(out) :: es(ii,kk) ! Saturation vapor pressure |
---|
265 | real(r8), intent(out) :: qs(ii,kk) ! Saturation specific humidity |
---|
266 | ! |
---|
267 | !---------------------------Local workspace----------------------------- |
---|
268 | ! |
---|
269 | real(r8) omeps ! 1 - 0.622 |
---|
270 | integer i, k ! Indices |
---|
271 | ! |
---|
272 | !----------------------------------------------------------------------- |
---|
273 | ! |
---|
274 | omeps = 1.0_r8 - epsqs |
---|
275 | do k=kstart,kend |
---|
276 | do i=1,ilen |
---|
277 | es(i,k) = estblf(t(i,k)) |
---|
278 | ! |
---|
279 | ! Saturation specific humidity |
---|
280 | ! |
---|
281 | qs(i,k) = epsqs*es(i,k)/(p(i,k) - omeps*es(i,k)) |
---|
282 | ! |
---|
283 | ! The following check is to avoid the generation of negative values |
---|
284 | ! that can occur in the upper stratosphere and mesosphere |
---|
285 | ! |
---|
286 | qs(i,k) = min(1.0_r8,qs(i,k)) |
---|
287 | ! |
---|
288 | if (qs(i,k) < 0.0_r8) then |
---|
289 | qs(i,k) = 1.0_r8 |
---|
290 | es(i,k) = p(i,k) |
---|
291 | end if |
---|
292 | end do |
---|
293 | end do |
---|
294 | ! |
---|
295 | return |
---|
296 | end subroutine aqsat |
---|
297 | |
---|
298 | #ifndef WRF_PORT |
---|
299 | !++xl |
---|
300 | subroutine aqsat_water(t ,p ,es ,qs ,ii , & |
---|
301 | ilen ,kk ,kstart ,kend ) |
---|
302 | !----------------------------------------------------------------------- |
---|
303 | ! |
---|
304 | ! Purpose: |
---|
305 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
306 | ! precomputed table, calculate and return saturation specific humidity |
---|
307 | ! (g/g),for input arrays of temperature and pressure (dimensioned ii,kk) |
---|
308 | ! This routine is useful for evaluating only a selected region in the |
---|
309 | ! vertical. |
---|
310 | ! |
---|
311 | ! Method: |
---|
312 | ! <Describe the algorithm(s) used in the routine.> |
---|
313 | ! <Also include any applicable external references.> |
---|
314 | ! |
---|
315 | ! Author: J. Hack |
---|
316 | ! |
---|
317 | !------------------------------Arguments-------------------------------- |
---|
318 | ! |
---|
319 | ! Input arguments |
---|
320 | ! |
---|
321 | integer, intent(in) :: ii ! I dimension of arrays t, p, es, qs |
---|
322 | integer, intent(in) :: kk ! K dimension of arrays t, p, es, qs |
---|
323 | integer, intent(in) :: ilen ! Length of vectors in I direction which |
---|
324 | integer, intent(in) :: kstart ! Starting location in K direction |
---|
325 | integer, intent(in) :: kend ! Ending location in K direction |
---|
326 | real(r8), intent(in) :: t(ii,kk) ! Temperature |
---|
327 | real(r8), intent(in) :: p(ii,kk) ! Pressure |
---|
328 | ! |
---|
329 | ! Output arguments |
---|
330 | ! |
---|
331 | real(r8), intent(out) :: es(ii,kk) ! Saturation vapor pressure |
---|
332 | real(r8), intent(out) :: qs(ii,kk) ! Saturation specific humidity |
---|
333 | ! |
---|
334 | !---------------------------Local workspace----------------------------- |
---|
335 | ! |
---|
336 | real(r8) omeps ! 1 - 0.622 |
---|
337 | integer i, k ! Indices |
---|
338 | ! |
---|
339 | !----------------------------------------------------------------------- |
---|
340 | ! |
---|
341 | omeps = 1.0_r8 - epsqs |
---|
342 | do k=kstart,kend |
---|
343 | do i=1,ilen |
---|
344 | ! es(i,k) = estblf(t(i,k)) |
---|
345 | es(i,k) = polysvp(t(i,k),0) |
---|
346 | ! |
---|
347 | ! Saturation specific humidity |
---|
348 | ! |
---|
349 | qs(i,k) = epsqs*es(i,k)/(p(i,k) - omeps*es(i,k)) |
---|
350 | ! |
---|
351 | ! The following check is to avoid the generation of negative values |
---|
352 | ! that can occur in the upper stratosphere and mesosphere |
---|
353 | ! |
---|
354 | qs(i,k) = min(1.0_r8,qs(i,k)) |
---|
355 | ! |
---|
356 | if (qs(i,k) < 0.0_r8) then |
---|
357 | qs(i,k) = 1.0_r8 |
---|
358 | es(i,k) = p(i,k) |
---|
359 | end if |
---|
360 | end do |
---|
361 | end do |
---|
362 | ! |
---|
363 | return |
---|
364 | end subroutine aqsat_water |
---|
365 | !--xl |
---|
366 | |
---|
367 | |
---|
368 | subroutine aqsatd(t ,p ,es ,qs ,gam , & |
---|
369 | ii ,ilen ,kk ,kstart ,kend ) |
---|
370 | !----------------------------------------------------------------------- |
---|
371 | ! |
---|
372 | ! Purpose: |
---|
373 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
374 | ! precomputed table, calculate and return saturation specific humidity |
---|
375 | ! (g/g). |
---|
376 | ! |
---|
377 | ! Method: |
---|
378 | ! Differs from aqsat by also calculating and returning |
---|
379 | ! gamma (l/cp)*(d(qsat)/dT) |
---|
380 | ! Input arrays temperature and pressure (dimensioned ii,kk). |
---|
381 | ! |
---|
382 | ! Author: J. Hack |
---|
383 | ! |
---|
384 | !------------------------------Arguments-------------------------------- |
---|
385 | ! |
---|
386 | ! Input arguments |
---|
387 | ! |
---|
388 | integer, intent(in) :: ii ! I dimension of arrays t, p, es, qs |
---|
389 | integer, intent(in) :: ilen ! Vector length in I direction |
---|
390 | integer, intent(in) :: kk ! K dimension of arrays t, p, es, qs |
---|
391 | integer, intent(in) :: kstart ! Starting location in K direction |
---|
392 | integer, intent(in) :: kend ! Ending location in K direction |
---|
393 | |
---|
394 | real(r8), intent(in) :: t(ii,kk) ! Temperature |
---|
395 | real(r8), intent(in) :: p(ii,kk) ! Pressure |
---|
396 | |
---|
397 | ! |
---|
398 | ! Output arguments |
---|
399 | ! |
---|
400 | real(r8), intent(out) :: es(ii,kk) ! Saturation vapor pressure |
---|
401 | real(r8), intent(out) :: qs(ii,kk) ! Saturation specific humidity |
---|
402 | real(r8), intent(out) :: gam(ii,kk) ! (l/cp)*(d(qs)/dt) |
---|
403 | ! |
---|
404 | !---------------------------Local workspace----------------------------- |
---|
405 | ! |
---|
406 | logical lflg ! True if in temperature transition region |
---|
407 | integer i ! i index for vector calculations |
---|
408 | integer k ! k index |
---|
409 | real(r8) omeps ! 1. - 0.622 |
---|
410 | real(r8) trinv ! Reciprocal of ttrice (transition range) |
---|
411 | real(r8) tc ! Temperature (in degrees C) |
---|
412 | real(r8) weight ! Weight for es transition from water to ice |
---|
413 | real(r8) hltalt ! Appropriately modified hlat for T derivatives |
---|
414 | real(r8) hlatsb ! hlat weighted in transition region |
---|
415 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
416 | real(r8) tterm ! Account for d(es)/dT in transition region |
---|
417 | real(r8) desdt ! d(es)/dT |
---|
418 | ! |
---|
419 | !----------------------------------------------------------------------- |
---|
420 | ! |
---|
421 | omeps = 1.0_r8 - epsqs |
---|
422 | do k=kstart,kend |
---|
423 | do i=1,ilen |
---|
424 | es(i,k) = estblf(t(i,k)) |
---|
425 | ! |
---|
426 | ! Saturation specific humidity |
---|
427 | ! |
---|
428 | qs(i,k) = epsqs*es(i,k)/(p(i,k) - omeps*es(i,k)) |
---|
429 | ! |
---|
430 | ! The following check is to avoid the generation of negative qs |
---|
431 | ! values which can occur in the upper stratosphere and mesosphere |
---|
432 | ! |
---|
433 | qs(i,k) = min(1.0_r8,qs(i,k)) |
---|
434 | ! |
---|
435 | if (qs(i,k) < 0.0_r8) then |
---|
436 | qs(i,k) = 1.0_r8 |
---|
437 | es(i,k) = p(i,k) |
---|
438 | end if |
---|
439 | end do |
---|
440 | end do |
---|
441 | ! |
---|
442 | ! "generalized" analytic expression for t derivative of es |
---|
443 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
444 | ! |
---|
445 | trinv = 0.0_r8 |
---|
446 | if ((.not. icephs) .or. (ttrice.eq.0.0_r8)) go to 10 |
---|
447 | trinv = 1.0_r8/ttrice |
---|
448 | ! |
---|
449 | do k=kstart,kend |
---|
450 | do i=1,ilen |
---|
451 | ! |
---|
452 | ! Weighting of hlat accounts for transition from water to ice |
---|
453 | ! polynomial expression approximates difference between es over |
---|
454 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
455 | ! -40): required for accurate estimate of es derivative in transition |
---|
456 | ! range from ice to water also accounting for change of hlatv with t |
---|
457 | ! above freezing where constant slope is given by -2369 j/(kg c) =cpv - cw |
---|
458 | ! |
---|
459 | tc = t(i,k) - tmelt |
---|
460 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
461 | weight = min(-tc*trinv,1.0_r8) |
---|
462 | hlatsb = hlatv + weight*hlatf |
---|
463 | hlatvp = hlatv - 2369.0_r8*tc |
---|
464 | if (t(i,k) < tmelt) then |
---|
465 | hltalt = hlatsb |
---|
466 | else |
---|
467 | hltalt = hlatvp |
---|
468 | end if |
---|
469 | if (lflg) then |
---|
470 | tterm = pcf(1) + tc*(pcf(2) + tc*(pcf(3) + tc*(pcf(4) + tc*pcf(5)))) |
---|
471 | else |
---|
472 | tterm = 0.0_r8 |
---|
473 | end if |
---|
474 | desdt = hltalt*es(i,k)/(rgasv*t(i,k)*t(i,k)) + tterm*trinv |
---|
475 | gam(i,k) = hltalt*qs(i,k)*p(i,k)*desdt/(cp*es(i,k)*(p(i,k) - omeps*es(i,k))) |
---|
476 | if (qs(i,k) == 1.0_r8) gam(i,k) = 0.0_r8 |
---|
477 | end do |
---|
478 | end do |
---|
479 | ! |
---|
480 | go to 20 |
---|
481 | ! |
---|
482 | ! No icephs or water to ice transition |
---|
483 | ! |
---|
484 | 10 do k=kstart,kend |
---|
485 | do i=1,ilen |
---|
486 | ! |
---|
487 | ! Account for change of hlatv with t above freezing where |
---|
488 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
489 | ! |
---|
490 | hlatvp = hlatv - 2369.0_r8*(t(i,k)-tmelt) |
---|
491 | if (icephs) then |
---|
492 | hlatsb = hlatv + hlatf |
---|
493 | else |
---|
494 | hlatsb = hlatv |
---|
495 | end if |
---|
496 | if (t(i,k) < tmelt) then |
---|
497 | hltalt = hlatsb |
---|
498 | else |
---|
499 | hltalt = hlatvp |
---|
500 | end if |
---|
501 | desdt = hltalt*es(i,k)/(rgasv*t(i,k)*t(i,k)) |
---|
502 | gam(i,k) = hltalt*qs(i,k)*p(i,k)*desdt/(cp*es(i,k)*(p(i,k) - omeps*es(i,k))) |
---|
503 | if (qs(i,k) == 1.0_r8) gam(i,k) = 0.0_r8 |
---|
504 | end do |
---|
505 | end do |
---|
506 | ! |
---|
507 | 20 return |
---|
508 | end subroutine aqsatd |
---|
509 | #endif |
---|
510 | |
---|
511 | subroutine vqsatd(t ,p ,es ,qs ,gam , & |
---|
512 | len ) |
---|
513 | !----------------------------------------------------------------------- |
---|
514 | ! |
---|
515 | ! Purpose: |
---|
516 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
517 | ! precomputed table, calculate and return saturation specific humidity |
---|
518 | ! (g/g), and calculate and return gamma (l/cp)*(d(qsat)/dT). The same |
---|
519 | ! function as qsatd, but operates on vectors of temperature and pressure |
---|
520 | ! |
---|
521 | ! Method: |
---|
522 | ! |
---|
523 | ! Author: J. Hack |
---|
524 | ! |
---|
525 | !------------------------------Arguments-------------------------------- |
---|
526 | ! |
---|
527 | ! Input arguments |
---|
528 | ! |
---|
529 | integer, intent(in) :: len ! vector length |
---|
530 | real(r8), intent(in) :: t(len) ! temperature |
---|
531 | real(r8), intent(in) :: p(len) ! pressure |
---|
532 | ! |
---|
533 | ! Output arguments |
---|
534 | ! |
---|
535 | real(r8), intent(out) :: es(len) ! saturation vapor pressure |
---|
536 | real(r8), intent(out) :: qs(len) ! saturation specific humidity |
---|
537 | real(r8), intent(out) :: gam(len) ! (l/cp)*(d(qs)/dt) |
---|
538 | ! |
---|
539 | !--------------------------Local Variables------------------------------ |
---|
540 | ! |
---|
541 | logical lflg ! true if in temperature transition region |
---|
542 | ! |
---|
543 | integer i ! index for vector calculations |
---|
544 | ! |
---|
545 | real(r8) omeps ! 1. - 0.622 |
---|
546 | real(r8) trinv ! reciprocal of ttrice (transition range) |
---|
547 | real(r8) tc ! temperature (in degrees C) |
---|
548 | real(r8) weight ! weight for es transition from water to ice |
---|
549 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
550 | ! |
---|
551 | real(r8) hlatsb ! hlat weighted in transition region |
---|
552 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
553 | real(r8) tterm ! account for d(es)/dT in transition region |
---|
554 | real(r8) desdt ! d(es)/dT |
---|
555 | ! |
---|
556 | !----------------------------------------------------------------------- |
---|
557 | ! |
---|
558 | omeps = 1.0_r8 - epsqs |
---|
559 | do i=1,len |
---|
560 | es(i) = estblf(t(i)) |
---|
561 | ! |
---|
562 | ! Saturation specific humidity |
---|
563 | ! |
---|
564 | qs(i) = epsqs*es(i)/(p(i) - omeps*es(i)) |
---|
565 | ! |
---|
566 | ! The following check is to avoid the generation of negative |
---|
567 | ! values that can occur in the upper stratosphere and mesosphere |
---|
568 | ! |
---|
569 | qs(i) = min(1.0_r8,qs(i)) |
---|
570 | ! |
---|
571 | if (qs(i) < 0.0_r8) then |
---|
572 | qs(i) = 1.0_r8 |
---|
573 | es(i) = p(i) |
---|
574 | end if |
---|
575 | end do |
---|
576 | ! |
---|
577 | ! "generalized" analytic expression for t derivative of es |
---|
578 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
579 | ! |
---|
580 | trinv = 0.0_r8 |
---|
581 | if ((.not. icephs) .or. (ttrice.eq.0.0_r8)) go to 10 |
---|
582 | trinv = 1.0_r8/ttrice |
---|
583 | do i=1,len |
---|
584 | ! |
---|
585 | ! Weighting of hlat accounts for transition from water to ice |
---|
586 | ! polynomial expression approximates difference between es over |
---|
587 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
588 | ! -40): required for accurate estimate of es derivative in transition |
---|
589 | ! range from ice to water also accounting for change of hlatv with t |
---|
590 | ! above freezing where const slope is given by -2369 j/(kg c) = cpv - cw |
---|
591 | ! |
---|
592 | tc = t(i) - tmelt |
---|
593 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
594 | weight = min(-tc*trinv,1.0_r8) |
---|
595 | hlatsb = hlatv + weight*hlatf |
---|
596 | hlatvp = hlatv - 2369.0_r8*tc |
---|
597 | if (t(i) < tmelt) then |
---|
598 | hltalt = hlatsb |
---|
599 | else |
---|
600 | hltalt = hlatvp |
---|
601 | end if |
---|
602 | if (lflg) then |
---|
603 | tterm = pcf(1) + tc*(pcf(2) + tc*(pcf(3) + tc*(pcf(4) + tc*pcf(5)))) |
---|
604 | else |
---|
605 | tterm = 0.0_r8 |
---|
606 | end if |
---|
607 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) + tterm*trinv |
---|
608 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
609 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
610 | end do |
---|
611 | return |
---|
612 | ! |
---|
613 | ! No icephs or water to ice transition |
---|
614 | ! |
---|
615 | 10 do i=1,len |
---|
616 | ! |
---|
617 | ! Account for change of hlatv with t above freezing where |
---|
618 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
619 | ! |
---|
620 | hlatvp = hlatv - 2369.0_r8*(t(i)-tmelt) |
---|
621 | if (icephs) then |
---|
622 | hlatsb = hlatv + hlatf |
---|
623 | else |
---|
624 | hlatsb = hlatv |
---|
625 | end if |
---|
626 | if (t(i) < tmelt) then |
---|
627 | hltalt = hlatsb |
---|
628 | else |
---|
629 | hltalt = hlatvp |
---|
630 | end if |
---|
631 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) |
---|
632 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
633 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
634 | end do |
---|
635 | ! |
---|
636 | return |
---|
637 | ! |
---|
638 | end subroutine vqsatd |
---|
639 | |
---|
640 | #ifndef WRF_PORT |
---|
641 | !++xl |
---|
642 | subroutine vqsatd_water(t ,p ,es ,qs ,gam , & |
---|
643 | len ) |
---|
644 | |
---|
645 | !------------------------------Arguments-------------------------------- |
---|
646 | ! |
---|
647 | ! Input arguments |
---|
648 | ! |
---|
649 | integer, intent(in) :: len ! vector length |
---|
650 | real(r8), intent(in) :: t(len) ! temperature |
---|
651 | real(r8), intent(in) :: p(len) ! pressure |
---|
652 | |
---|
653 | ! |
---|
654 | ! Output arguments |
---|
655 | ! |
---|
656 | real(r8), intent(out) :: es(len) ! saturation vapor pressure |
---|
657 | real(r8), intent(out) :: qs(len) ! saturation specific humidity |
---|
658 | real(r8), intent(out) :: gam(len) ! (l/cp)*(d(qs)/dt) |
---|
659 | |
---|
660 | ! |
---|
661 | !--------------------------Local Variables------------------------------ |
---|
662 | ! |
---|
663 | ! |
---|
664 | integer i ! index for vector calculations |
---|
665 | ! |
---|
666 | real(r8) omeps ! 1. - 0.622 |
---|
667 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
668 | ! |
---|
669 | real(r8) hlatsb ! hlat weighted in transition region |
---|
670 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
671 | real(r8) desdt ! d(es)/dT |
---|
672 | ! |
---|
673 | !----------------------------------------------------------------------- |
---|
674 | ! |
---|
675 | omeps = 1.0_r8 - epsqs |
---|
676 | do i=1,len |
---|
677 | es(i) = polysvp(t(i),0) |
---|
678 | ! |
---|
679 | ! Saturation specific humidity |
---|
680 | ! |
---|
681 | qs(i) = epsqs*es(i)/(p(i) - omeps*es(i)) |
---|
682 | ! |
---|
683 | ! The following check is to avoid the generation of negative |
---|
684 | ! values that can occur in the upper stratosphere and mesosphere |
---|
685 | ! |
---|
686 | qs(i) = min(1.0_r8,qs(i)) |
---|
687 | ! |
---|
688 | if (qs(i) < 0.0_r8) then |
---|
689 | qs(i) = 1.0_r8 |
---|
690 | es(i) = p(i) |
---|
691 | end if |
---|
692 | end do |
---|
693 | ! |
---|
694 | ! No icephs or water to ice transition |
---|
695 | ! |
---|
696 | do i=1,len |
---|
697 | ! |
---|
698 | ! Account for change of hlatv with t above freezing where |
---|
699 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
700 | ! |
---|
701 | hlatvp = hlatv - 2369.0_r8*(t(i)-tmelt) |
---|
702 | hlatsb = hlatv |
---|
703 | if (t(i) < tmelt) then |
---|
704 | hltalt = hlatsb |
---|
705 | else |
---|
706 | hltalt = hlatvp |
---|
707 | end if |
---|
708 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) |
---|
709 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
710 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
711 | end do |
---|
712 | ! |
---|
713 | return |
---|
714 | ! |
---|
715 | end subroutine vqsatd_water |
---|
716 | |
---|
717 | function polysvp (T,type) |
---|
718 | ! Compute saturation vapor pressure by using |
---|
719 | ! function from Goff and Gatch (1946) |
---|
720 | |
---|
721 | ! Polysvp returned in units of pa. |
---|
722 | ! T is input in units of K. |
---|
723 | ! type refers to saturation with respect to liquid (0) or ice (1) |
---|
724 | |
---|
725 | real(r8) dum |
---|
726 | |
---|
727 | real(r8) T,polysvp |
---|
728 | |
---|
729 | integer type |
---|
730 | |
---|
731 | ! ice |
---|
732 | |
---|
733 | if (type.eq.1) then |
---|
734 | |
---|
735 | ! Goff Gatch equation (good down to -100 C) |
---|
736 | |
---|
737 | polysvp = 10._r8**(-9.09718_r8*(273.16_r8/t-1._r8)-3.56654_r8* & |
---|
738 | log10(273.16_r8/t)+0.876793_r8*(1._r8-t/273.16_r8)+ & |
---|
739 | log10(6.1071_r8))*100._r8 |
---|
740 | |
---|
741 | end if |
---|
742 | |
---|
743 | ! Goff Gatch equation, uncertain below -70 C |
---|
744 | |
---|
745 | if (type.eq.0) then |
---|
746 | polysvp = 10._r8**(-7.90298_r8*(373.16_r8/t-1._r8)+ & |
---|
747 | 5.02808_r8*log10(373.16_r8/t)- & |
---|
748 | 1.3816e-7_r8*(10._r8**(11.344_r8*(1._r8-t/373.16_r8))-1._r8)+ & |
---|
749 | 8.1328e-3_r8*(10._r8**(-3.49149_r8*(373.16_r8/t-1._r8))-1._r8)+ & |
---|
750 | log10(1013.246_r8))*100._r8 |
---|
751 | end if |
---|
752 | |
---|
753 | |
---|
754 | end function polysvp |
---|
755 | !--xl |
---|
756 | #endif |
---|
757 | |
---|
758 | integer function fqsatd(t ,p ,es ,qs ,gam , len ) |
---|
759 | !----------------------------------------------------------------------- |
---|
760 | ! Purpose: |
---|
761 | ! This is merely a function interface vqsatd. |
---|
762 | !------------------------------Arguments-------------------------------- |
---|
763 | ! Input arguments |
---|
764 | integer, intent(in) :: len ! vector length |
---|
765 | real(r8), intent(in) :: t(len) ! temperature |
---|
766 | real(r8), intent(in) :: p(len) ! pressure |
---|
767 | ! Output arguments |
---|
768 | real(r8), intent(out) :: es(len) ! saturation vapor pressure |
---|
769 | real(r8), intent(out) :: qs(len) ! saturation specific humidity |
---|
770 | real(r8), intent(out) :: gam(len) ! (l/cp)*(d(qs)/dt) |
---|
771 | ! Call vqsatd |
---|
772 | call vqsatd(t ,p ,es ,qs ,gam , len ) |
---|
773 | fqsatd = 1 |
---|
774 | return |
---|
775 | end function fqsatd |
---|
776 | |
---|
777 | #ifndef WRF_PORT |
---|
778 | real(r8) function qsat_water(t,p) |
---|
779 | ! saturation mixing ratio w/respect to liquid water |
---|
780 | real(r8) t ! temperature |
---|
781 | real(r8) p ! pressure (Pa) |
---|
782 | real(r8) es ! saturation vapor pressure (Pa) |
---|
783 | real(r8) ps, ts, e1, e2, f1, f2, f3, f4, f5, f |
---|
784 | ! real(r8) t0inv ! 1/273. |
---|
785 | ! data t0inv/0.003663/ |
---|
786 | ! save t0inv |
---|
787 | ! es = 611.*exp(hlatv/rgasv*(t0inv-1./t)) |
---|
788 | |
---|
789 | ps = 1013.246_r8 |
---|
790 | ts = 373.16_r8 |
---|
791 | e1 = 11.344_r8*(1.0_r8 - t/ts) |
---|
792 | e2 = -3.49149_r8*(ts/t - 1.0_r8) |
---|
793 | f1 = -7.90298_r8*(ts/t - 1.0_r8) |
---|
794 | f2 = 5.02808_r8*log10(ts/t) |
---|
795 | f3 = -1.3816_r8*(10.0_r8**e1 - 1.0_r8)/10000000.0_r8 |
---|
796 | f4 = 8.1328_r8*(10.0_r8**e2 - 1.0_r8)/1000.0_r8 |
---|
797 | f5 = log10(ps) |
---|
798 | f = f1 + f2 + f3 + f4 + f5 |
---|
799 | es = (10.0_r8**f)*100.0_r8 |
---|
800 | |
---|
801 | qsat_water = epsqs*es/(p-(1.-epsqs)*es) ! saturation w/respect to liquid only |
---|
802 | if(qsat_water < 0.) qsat_water = 1. |
---|
803 | |
---|
804 | end function qsat_water |
---|
805 | |
---|
806 | subroutine vqsat_water(t,p,qsat_water,len) |
---|
807 | ! saturation mixing ratio w/respect to liquid water |
---|
808 | integer, intent(in) :: len |
---|
809 | real(r8) t(len) ! temperature |
---|
810 | real(r8) p(len) ! pressure (Pa) |
---|
811 | real(r8) qsat_water(len) |
---|
812 | real(r8) es ! saturation vapor pressure (Pa) |
---|
813 | real(r8), parameter :: t0inv = 1._r8/273._r8 |
---|
814 | real(r8) coef |
---|
815 | integer :: i |
---|
816 | |
---|
817 | coef = hlatv/rgasv |
---|
818 | do i=1,len |
---|
819 | es = 611._r8*exp(coef*(t0inv-1./t(i))) |
---|
820 | qsat_water(i) = epsqs*es/(p(i)-(1.-epsqs)*es) ! saturation w/respect to liquid only |
---|
821 | if(qsat_water(i) < 0.) qsat_water(i) = 1. |
---|
822 | enddo |
---|
823 | |
---|
824 | return |
---|
825 | |
---|
826 | end subroutine vqsat_water |
---|
827 | |
---|
828 | real(r8) function qsat_ice(t,p) |
---|
829 | ! saturation mixing ratio w/respect to ice |
---|
830 | real(r8) t ! temperature |
---|
831 | real(r8) p ! pressure (Pa) |
---|
832 | real(r8) es ! saturation vapor pressure (Pa) |
---|
833 | real(r8), parameter :: t0inv = 1._r8/273._r8 |
---|
834 | es = 611.*exp((hlatv+hlatf)/rgasv*(t0inv-1./t)) |
---|
835 | qsat_ice = epsqs*es/(p-(1.-epsqs)*es) ! saturation w/respect to liquid only |
---|
836 | if(qsat_ice < 0.) qsat_ice = 1. |
---|
837 | |
---|
838 | end function qsat_ice |
---|
839 | |
---|
840 | subroutine vqsat_ice(t,p,qsat_ice,len) |
---|
841 | ! saturation mixing ratio w/respect to liquid water |
---|
842 | integer,intent(in) :: len |
---|
843 | real(r8) t(len) ! temperature |
---|
844 | real(r8) p(len) ! pressure (Pa) |
---|
845 | real(r8) qsat_ice(len) |
---|
846 | real(r8) es ! saturation vapor pressure (Pa) |
---|
847 | real(r8), parameter :: t0inv = 1._r8/273._r8 |
---|
848 | real(r8) coef |
---|
849 | integer :: i |
---|
850 | |
---|
851 | coef = (hlatv+hlatf)/rgasv |
---|
852 | do i=1,len |
---|
853 | es = 611.*exp(coef*(t0inv-1./t(i))) |
---|
854 | qsat_ice(i) = epsqs*es/(p(i)-(1.-epsqs)*es) ! saturation w/respect to liquid only |
---|
855 | if(qsat_ice(i) < 0.) qsat_ice(i) = 1. |
---|
856 | enddo |
---|
857 | |
---|
858 | return |
---|
859 | |
---|
860 | end subroutine vqsat_ice |
---|
861 | |
---|
862 | ! Sungsu |
---|
863 | ! Below two subroutines (vqsatd2_water,vqsatd2_water_single) are by Sungsu |
---|
864 | ! Replace 'gam -> dqsdt' |
---|
865 | ! Sungsu |
---|
866 | |
---|
867 | subroutine vqsatd2_water(t ,p ,es ,qs ,dqsdt , & |
---|
868 | len ) |
---|
869 | |
---|
870 | !------------------------------Arguments-------------------------------- |
---|
871 | ! |
---|
872 | ! Input arguments |
---|
873 | ! |
---|
874 | integer, intent(in) :: len ! vector length |
---|
875 | real(r8), intent(in) :: t(len) ! temperature |
---|
876 | real(r8), intent(in) :: p(len) ! pressure |
---|
877 | |
---|
878 | ! |
---|
879 | ! Output arguments |
---|
880 | ! |
---|
881 | real(r8), intent(out) :: es(len) ! saturation vapor pressure |
---|
882 | real(r8), intent(out) :: qs(len) ! saturation specific humidity |
---|
883 | ! real(r8), intent(out) :: gam(len) ! (l/cp)*(d(qs)/dt) |
---|
884 | ! Sungsu |
---|
885 | real(r8), intent(out) :: dqsdt(len) ! (d(qs)/dt) |
---|
886 | ! End by Sungsu |
---|
887 | |
---|
888 | ! |
---|
889 | !--------------------------Local Variables------------------------------ |
---|
890 | ! |
---|
891 | ! |
---|
892 | integer i ! index for vector calculations |
---|
893 | ! |
---|
894 | real(r8) omeps ! 1. - 0.622 |
---|
895 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
896 | ! |
---|
897 | real(r8) hlatsb ! hlat weighted in transition region |
---|
898 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
899 | real(r8) desdt ! d(es)/dT |
---|
900 | |
---|
901 | ! Sungsu |
---|
902 | real(r8) gam(len) ! (l/cp)*(d(qs)/dt) |
---|
903 | ! End by Sungsu |
---|
904 | |
---|
905 | ! |
---|
906 | !----------------------------------------------------------------------- |
---|
907 | ! |
---|
908 | omeps = 1.0_r8 - epsqs |
---|
909 | do i=1,len |
---|
910 | es(i) = polysvp(t(i),0) |
---|
911 | ! |
---|
912 | ! Saturation specific humidity |
---|
913 | ! |
---|
914 | qs(i) = epsqs*es(i)/(p(i) - omeps*es(i)) |
---|
915 | ! |
---|
916 | ! The following check is to avoid the generation of negative |
---|
917 | ! values that can occur in the upper stratosphere and mesosphere |
---|
918 | ! |
---|
919 | qs(i) = min(1.0_r8,qs(i)) |
---|
920 | ! |
---|
921 | if (qs(i) < 0.0_r8) then |
---|
922 | qs(i) = 1.0_r8 |
---|
923 | es(i) = p(i) |
---|
924 | end if |
---|
925 | end do |
---|
926 | ! |
---|
927 | ! No icephs or water to ice transition |
---|
928 | ! |
---|
929 | do i=1,len |
---|
930 | ! |
---|
931 | ! Account for change of hlatv with t above freezing where |
---|
932 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
933 | ! |
---|
934 | hlatvp = hlatv - 2369.0_r8*(t(i)-tmelt) |
---|
935 | hlatsb = hlatv |
---|
936 | if (t(i) < tmelt) then |
---|
937 | hltalt = hlatsb |
---|
938 | else |
---|
939 | hltalt = hlatvp |
---|
940 | end if |
---|
941 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) |
---|
942 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
943 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
944 | ! Sungsu |
---|
945 | dqsdt(i) = (cp/hltalt)*gam(i) |
---|
946 | ! End by Sungsu |
---|
947 | end do |
---|
948 | ! |
---|
949 | return |
---|
950 | ! |
---|
951 | end subroutine vqsatd2_water |
---|
952 | |
---|
953 | subroutine vqsatd2_water_single(t ,p ,es ,qs ,dqsdt) |
---|
954 | |
---|
955 | !------------------------------Arguments-------------------------------- |
---|
956 | ! |
---|
957 | ! Input arguments |
---|
958 | ! |
---|
959 | |
---|
960 | real(r8), intent(in) :: t ! temperature |
---|
961 | real(r8), intent(in) :: p ! pressure |
---|
962 | |
---|
963 | ! |
---|
964 | ! Output arguments |
---|
965 | ! |
---|
966 | real(r8), intent(out) :: es ! saturation vapor pressure |
---|
967 | real(r8), intent(out) :: qs ! saturation specific humidity |
---|
968 | ! real(r8), intent(out) :: gam ! (l/cp)*(d(qs)/dt) |
---|
969 | ! Sungsu |
---|
970 | real(r8), intent(out) :: dqsdt ! (d(qs)/dt) |
---|
971 | ! End by Sungsu |
---|
972 | |
---|
973 | ! |
---|
974 | !--------------------------Local Variables------------------------------ |
---|
975 | ! |
---|
976 | ! |
---|
977 | integer i ! index for vector calculations |
---|
978 | ! |
---|
979 | real(r8) omeps ! 1. - 0.622 |
---|
980 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
981 | ! |
---|
982 | real(r8) hlatsb ! hlat weighted in transition region |
---|
983 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
984 | real(r8) desdt ! d(es)/dT |
---|
985 | |
---|
986 | ! Sungsu |
---|
987 | real(r8) gam ! (l/cp)*(d(qs)/dt) |
---|
988 | ! End by Sungsu |
---|
989 | |
---|
990 | ! |
---|
991 | !----------------------------------------------------------------------- |
---|
992 | ! |
---|
993 | omeps = 1.0_r8 - epsqs |
---|
994 | ! do i=1,len |
---|
995 | es = polysvp(t,0) |
---|
996 | ! |
---|
997 | ! Saturation specific humidity |
---|
998 | ! |
---|
999 | qs = epsqs*es/(p - omeps*es) |
---|
1000 | ! |
---|
1001 | ! The following check is to avoid the generation of negative |
---|
1002 | ! values that can occur in the upper stratosphere and mesosphere |
---|
1003 | ! |
---|
1004 | qs = min(1.0_r8,qs) |
---|
1005 | ! |
---|
1006 | if (qs < 0.0_r8) then |
---|
1007 | qs = 1.0_r8 |
---|
1008 | es = p |
---|
1009 | end if |
---|
1010 | ! end do |
---|
1011 | ! |
---|
1012 | ! No icephs or water to ice transition |
---|
1013 | ! |
---|
1014 | ! do i=1,len |
---|
1015 | ! |
---|
1016 | ! Account for change of hlatv with t above freezing where |
---|
1017 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
1018 | ! |
---|
1019 | hlatvp = hlatv - 2369.0_r8*(t-tmelt) |
---|
1020 | hlatsb = hlatv |
---|
1021 | if (t < tmelt) then |
---|
1022 | hltalt = hlatsb |
---|
1023 | else |
---|
1024 | hltalt = hlatvp |
---|
1025 | end if |
---|
1026 | desdt = hltalt*es/(rgasv*t*t) |
---|
1027 | gam = hltalt*qs*p*desdt/(cp*es*(p - omeps*es)) |
---|
1028 | if (qs == 1.0_r8) gam = 0.0_r8 |
---|
1029 | ! Sungsu |
---|
1030 | dqsdt = (cp/hltalt)*gam |
---|
1031 | ! End by Sungsu |
---|
1032 | ! end do |
---|
1033 | ! |
---|
1034 | return |
---|
1035 | ! |
---|
1036 | end subroutine vqsatd2_water_single |
---|
1037 | |
---|
1038 | |
---|
1039 | subroutine vqsatd2(t ,p ,es ,qs ,dqsdt , & |
---|
1040 | len ) |
---|
1041 | !----------------------------------------------------------------------- |
---|
1042 | ! Sungsu : This is directly copied from 'vqsatd' but 'dqsdt' is output |
---|
1043 | ! instead of gam for use in Sungsu's equilibrium stratiform |
---|
1044 | ! macrophysics scheme. |
---|
1045 | ! |
---|
1046 | ! Purpose: |
---|
1047 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
1048 | ! precomputed table, calculate and return saturation specific humidity |
---|
1049 | ! (g/g), and calculate and return gamma (l/cp)*(d(qsat)/dT). The same |
---|
1050 | ! function as qsatd, but operates on vectors of temperature and pressure |
---|
1051 | ! |
---|
1052 | ! Method: |
---|
1053 | ! |
---|
1054 | ! Author: J. Hack |
---|
1055 | ! |
---|
1056 | !------------------------------Arguments-------------------------------- |
---|
1057 | ! |
---|
1058 | ! Input arguments |
---|
1059 | ! |
---|
1060 | integer, intent(in) :: len ! vector length |
---|
1061 | real(r8), intent(in) :: t(len) ! temperature |
---|
1062 | real(r8), intent(in) :: p(len) ! pressure |
---|
1063 | ! |
---|
1064 | ! Output arguments |
---|
1065 | ! |
---|
1066 | real(r8), intent(out) :: es(len) ! saturation vapor pressure |
---|
1067 | real(r8), intent(out) :: qs(len) ! saturation specific humidity |
---|
1068 | ! real(r8), intent(out) :: gam(len) ! (l/cp)*(d(qs)/dt) |
---|
1069 | ! Sungsu |
---|
1070 | real(r8), intent(out) :: dqsdt(len) ! (d(qs)/dt) |
---|
1071 | ! End by Sungsu |
---|
1072 | |
---|
1073 | ! |
---|
1074 | !--------------------------Local Variables------------------------------ |
---|
1075 | ! |
---|
1076 | logical lflg ! true if in temperature transition region |
---|
1077 | ! |
---|
1078 | integer i ! index for vector calculations |
---|
1079 | ! |
---|
1080 | real(r8) omeps ! 1. - 0.622 |
---|
1081 | real(r8) trinv ! reciprocal of ttrice (transition range) |
---|
1082 | real(r8) tc ! temperature (in degrees C) |
---|
1083 | real(r8) weight ! weight for es transition from water to ice |
---|
1084 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
1085 | ! |
---|
1086 | real(r8) hlatsb ! hlat weighted in transition region |
---|
1087 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
1088 | real(r8) tterm ! account for d(es)/dT in transition region |
---|
1089 | real(r8) desdt ! d(es)/dT |
---|
1090 | |
---|
1091 | ! Sungsu |
---|
1092 | real(r8) gam(len) ! (l/cp)*(d(qs)/dt) |
---|
1093 | ! End by Sungsu |
---|
1094 | ! |
---|
1095 | !----------------------------------------------------------------------- |
---|
1096 | ! |
---|
1097 | omeps = 1.0_r8 - epsqs |
---|
1098 | do i=1,len |
---|
1099 | es(i) = estblf(t(i)) |
---|
1100 | ! |
---|
1101 | ! Saturation specific humidity |
---|
1102 | ! |
---|
1103 | qs(i) = epsqs*es(i)/(p(i) - omeps*es(i)) |
---|
1104 | ! |
---|
1105 | ! The following check is to avoid the generation of negative |
---|
1106 | ! values that can occur in the upper stratosphere and mesosphere |
---|
1107 | ! |
---|
1108 | qs(i) = min(1.0_r8,qs(i)) |
---|
1109 | ! |
---|
1110 | if (qs(i) < 0.0_r8) then |
---|
1111 | qs(i) = 1.0_r8 |
---|
1112 | es(i) = p(i) |
---|
1113 | end if |
---|
1114 | end do |
---|
1115 | ! |
---|
1116 | ! "generalized" analytic expression for t derivative of es |
---|
1117 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
1118 | ! |
---|
1119 | trinv = 0.0_r8 |
---|
1120 | if ((.not. icephs) .or. (ttrice.eq.0.0_r8)) go to 10 |
---|
1121 | trinv = 1.0_r8/ttrice |
---|
1122 | do i=1,len |
---|
1123 | ! |
---|
1124 | ! Weighting of hlat accounts for transition from water to ice |
---|
1125 | ! polynomial expression approximates difference between es over |
---|
1126 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
1127 | ! -40): required for accurate estimate of es derivative in transition |
---|
1128 | ! range from ice to water also accounting for change of hlatv with t |
---|
1129 | ! above freezing where const slope is given by -2369 j/(kg c) = cpv - cw |
---|
1130 | ! |
---|
1131 | tc = t(i) - tmelt |
---|
1132 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
1133 | weight = min(-tc*trinv,1.0_r8) |
---|
1134 | hlatsb = hlatv + weight*hlatf |
---|
1135 | hlatvp = hlatv - 2369.0_r8*tc |
---|
1136 | if (t(i) < tmelt) then |
---|
1137 | hltalt = hlatsb |
---|
1138 | else |
---|
1139 | hltalt = hlatvp |
---|
1140 | end if |
---|
1141 | if (lflg) then |
---|
1142 | tterm = pcf(1) + tc*(pcf(2) + tc*(pcf(3) + tc*(pcf(4) + tc*pcf(5)))) |
---|
1143 | else |
---|
1144 | tterm = 0.0_r8 |
---|
1145 | end if |
---|
1146 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) + tterm*trinv |
---|
1147 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
1148 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
1149 | ! Sungsu |
---|
1150 | dqsdt(i) = (cp/hltalt)*gam(i) |
---|
1151 | ! End by Sungsu |
---|
1152 | end do |
---|
1153 | return |
---|
1154 | ! |
---|
1155 | ! No icephs or water to ice transition |
---|
1156 | ! |
---|
1157 | 10 do i=1,len |
---|
1158 | ! |
---|
1159 | ! Account for change of hlatv with t above freezing where |
---|
1160 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
1161 | ! |
---|
1162 | hlatvp = hlatv - 2369.0_r8*(t(i)-tmelt) |
---|
1163 | if (icephs) then |
---|
1164 | hlatsb = hlatv + hlatf |
---|
1165 | else |
---|
1166 | hlatsb = hlatv |
---|
1167 | end if |
---|
1168 | if (t(i) < tmelt) then |
---|
1169 | hltalt = hlatsb |
---|
1170 | else |
---|
1171 | hltalt = hlatvp |
---|
1172 | end if |
---|
1173 | desdt = hltalt*es(i)/(rgasv*t(i)*t(i)) |
---|
1174 | gam(i) = hltalt*qs(i)*p(i)*desdt/(cp*es(i)*(p(i) - omeps*es(i))) |
---|
1175 | if (qs(i) == 1.0_r8) gam(i) = 0.0_r8 |
---|
1176 | ! Sungsu |
---|
1177 | dqsdt(i) = (cp/hltalt)*gam(i) |
---|
1178 | ! End by Sungsu |
---|
1179 | end do |
---|
1180 | ! |
---|
1181 | return |
---|
1182 | ! |
---|
1183 | end subroutine vqsatd2 |
---|
1184 | |
---|
1185 | |
---|
1186 | ! Below routine is by Sungsu |
---|
1187 | |
---|
1188 | subroutine vqsatd2_single(t ,p ,es ,qs ,dqsdt) |
---|
1189 | !----------------------------------------------------------------------- |
---|
1190 | ! Sungsu : This is directly copied from 'vqsatd' but 'dqsdt' is output |
---|
1191 | ! instead of gam for use in Sungsu's equilibrium stratiform |
---|
1192 | ! macrophysics scheme. |
---|
1193 | ! |
---|
1194 | ! Purpose: |
---|
1195 | ! Utility procedure to look up and return saturation vapor pressure from |
---|
1196 | ! precomputed table, calculate and return saturation specific humidity |
---|
1197 | ! (g/g), and calculate and return gamma (l/cp)*(d(qsat)/dT). The same |
---|
1198 | ! function as qsatd, but operates on vectors of temperature and pressure |
---|
1199 | ! |
---|
1200 | ! Method: |
---|
1201 | ! |
---|
1202 | ! Author: J. Hack |
---|
1203 | ! |
---|
1204 | !------------------------------Arguments-------------------------------- |
---|
1205 | ! |
---|
1206 | ! Input arguments |
---|
1207 | ! |
---|
1208 | real(r8), intent(in) :: t ! temperature |
---|
1209 | real(r8), intent(in) :: p ! pressure |
---|
1210 | ! |
---|
1211 | ! Output arguments |
---|
1212 | ! |
---|
1213 | real(r8), intent(out) :: es ! saturation vapor pressure |
---|
1214 | real(r8), intent(out) :: qs ! saturation specific humidity |
---|
1215 | ! real(r8), intent(out) :: gam ! (l/cp)*(d(qs)/dt) |
---|
1216 | ! Sungsu |
---|
1217 | real(r8), intent(out) :: dqsdt ! (d(qs)/dt) |
---|
1218 | ! End by Sungsu |
---|
1219 | |
---|
1220 | ! |
---|
1221 | !--------------------------Local Variables------------------------------ |
---|
1222 | ! |
---|
1223 | logical lflg ! true if in temperature transition region |
---|
1224 | ! |
---|
1225 | ! integer i ! index for vector calculations |
---|
1226 | ! |
---|
1227 | real(r8) omeps ! 1. - 0.622 |
---|
1228 | real(r8) trinv ! reciprocal of ttrice (transition range) |
---|
1229 | real(r8) tc ! temperature (in degrees C) |
---|
1230 | real(r8) weight ! weight for es transition from water to ice |
---|
1231 | real(r8) hltalt ! appropriately modified hlat for T derivatives |
---|
1232 | ! |
---|
1233 | real(r8) hlatsb ! hlat weighted in transition region |
---|
1234 | real(r8) hlatvp ! hlat modified for t changes above freezing |
---|
1235 | real(r8) tterm ! account for d(es)/dT in transition region |
---|
1236 | real(r8) desdt ! d(es)/dT |
---|
1237 | |
---|
1238 | ! Sungsu |
---|
1239 | real(r8) gam ! (l/cp)*(d(qs)/dt) |
---|
1240 | ! End by Sungsu |
---|
1241 | ! |
---|
1242 | !----------------------------------------------------------------------- |
---|
1243 | ! |
---|
1244 | omeps = 1.0_r8 - epsqs |
---|
1245 | |
---|
1246 | ! do i=1,len |
---|
1247 | |
---|
1248 | es = estblf(t) |
---|
1249 | ! |
---|
1250 | ! Saturation specific humidity |
---|
1251 | ! |
---|
1252 | qs = epsqs*es/(p - omeps*es) |
---|
1253 | ! |
---|
1254 | ! The following check is to avoid the generation of negative |
---|
1255 | ! values that can occur in the upper stratosphere and mesosphere |
---|
1256 | ! |
---|
1257 | qs = min(1.0_r8,qs) |
---|
1258 | ! |
---|
1259 | if (qs < 0.0_r8) then |
---|
1260 | qs = 1.0_r8 |
---|
1261 | es = p |
---|
1262 | end if |
---|
1263 | |
---|
1264 | ! end do |
---|
1265 | ! |
---|
1266 | ! "generalized" analytic expression for t derivative of es |
---|
1267 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
1268 | ! |
---|
1269 | trinv = 0.0_r8 |
---|
1270 | if ((.not. icephs) .or. (ttrice.eq.0.0_r8)) go to 10 |
---|
1271 | trinv = 1.0_r8/ttrice |
---|
1272 | |
---|
1273 | ! do i=1,len |
---|
1274 | ! |
---|
1275 | ! Weighting of hlat accounts for transition from water to ice |
---|
1276 | ! polynomial expression approximates difference between es over |
---|
1277 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
1278 | ! -40): required for accurate estimate of es derivative in transition |
---|
1279 | ! range from ice to water also accounting for change of hlatv with t |
---|
1280 | ! above freezing where const slope is given by -2369 j/(kg c) = cpv - cw |
---|
1281 | ! |
---|
1282 | tc = t - tmelt |
---|
1283 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
1284 | weight = min(-tc*trinv,1.0_r8) |
---|
1285 | hlatsb = hlatv + weight*hlatf |
---|
1286 | hlatvp = hlatv - 2369.0_r8*tc |
---|
1287 | if (t < tmelt) then |
---|
1288 | hltalt = hlatsb |
---|
1289 | else |
---|
1290 | hltalt = hlatvp |
---|
1291 | end if |
---|
1292 | if (lflg) then |
---|
1293 | tterm = pcf(1) + tc*(pcf(2) + tc*(pcf(3) + tc*(pcf(4) + tc*pcf(5)))) |
---|
1294 | else |
---|
1295 | tterm = 0.0_r8 |
---|
1296 | end if |
---|
1297 | desdt = hltalt*es/(rgasv*t*t) + tterm*trinv |
---|
1298 | gam = hltalt*qs*p*desdt/(cp*es*(p - omeps*es)) |
---|
1299 | if (qs == 1.0_r8) gam = 0.0_r8 |
---|
1300 | ! Sungsu |
---|
1301 | dqsdt = (cp/hltalt)*gam |
---|
1302 | ! End by Sungsu |
---|
1303 | ! end do |
---|
1304 | return |
---|
1305 | ! |
---|
1306 | ! No icephs or water to ice transition |
---|
1307 | ! |
---|
1308 | |
---|
1309 | 10 continue |
---|
1310 | |
---|
1311 | !10 do i=1,len |
---|
1312 | ! |
---|
1313 | ! Account for change of hlatv with t above freezing where |
---|
1314 | ! constant slope is given by -2369 j/(kg c) = cpv - cw |
---|
1315 | ! |
---|
1316 | hlatvp = hlatv - 2369.0_r8*(t-tmelt) |
---|
1317 | if (icephs) then |
---|
1318 | hlatsb = hlatv + hlatf |
---|
1319 | else |
---|
1320 | hlatsb = hlatv |
---|
1321 | end if |
---|
1322 | if (t < tmelt) then |
---|
1323 | hltalt = hlatsb |
---|
1324 | else |
---|
1325 | hltalt = hlatvp |
---|
1326 | end if |
---|
1327 | desdt = hltalt*es/(rgasv*t*t) |
---|
1328 | gam = hltalt*qs*p*desdt/(cp*es*(p - omeps*es)) |
---|
1329 | if (qs == 1.0_r8) gam = 0.0_r8 |
---|
1330 | ! Sungsu |
---|
1331 | dqsdt = (cp/hltalt)*gam |
---|
1332 | ! End by Sungsu |
---|
1333 | |
---|
1334 | ! end do |
---|
1335 | ! |
---|
1336 | return |
---|
1337 | ! |
---|
1338 | end subroutine vqsatd2_single |
---|
1339 | #endif |
---|
1340 | |
---|
1341 | end module wv_saturation |
---|