Interoperable Physics Driver for NGGPS
funcphys.f
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 module funcphys
3 !$$$ Module Documentation Block
4 !
5 ! Module: funcphys API for basic thermodynamic physics
6 ! Author: Iredell Org: W/NX23 Date: 1999-03-01
7 !
8 ! Abstract: This module provides an Application Program Interface
9 ! for computing basic thermodynamic physics functions, in particular
10 ! (1) saturation vapor pressure as a function of temperature,
11 ! (2) dewpoint temperature as a function of vapor pressure,
12 ! (3) equivalent potential temperature as a function of temperature
13 ! and scaled pressure to the kappa power,
14 ! (4) temperature and specific humidity along a moist adiabat
15 ! as functions of equivalent potential temperature and
16 ! scaled pressure to the kappa power,
17 ! (5) scaled pressure to the kappa power as a function of pressure, and
18 ! (6) temperature at the lifting condensation level as a function
19 ! of temperature and dewpoint depression.
20 ! The entry points required to set up lookup tables start with a "g".
21 ! All the other entry points are functions starting with an "f" or
22 ! are subroutines starting with an "s". These other functions and
23 ! subroutines are elemental; that is, they return a scalar if they
24 ! are passed only scalars, but they return an array if they are passed
25 ! an array. These other functions and subroutines can be inlined, too.
26 !
27 ! Program History Log:
28 ! 1999-03-01 Mark Iredell
29 ! 1999-10-15 Mark Iredell SI unit for pressure (Pascals)
30 ! 2001-02-26 Mark Iredell Ice phase changes of Hong and Moorthi
31 !
32 ! Public Variables:
33 ! krealfp Integer parameter kind or length of reals (=kind_phys)
34 !
35 ! Public Subprograms:
36 ! gpvsl Compute saturation vapor pressure over liquid table
37 !
38 ! fpvsl Elementally compute saturation vapor pressure over liquid
39 ! function result Real(krealfp) saturation vapor pressure in Pascals
40 ! t Real(krealfp) temperature in Kelvin
41 !
42 ! fpvslq Elementally compute saturation vapor pressure over liquid
43 ! function result Real(krealfp) saturation vapor pressure in Pascals
44 ! t Real(krealfp) temperature in Kelvin
45 !
46 ! fpvslx Elementally compute saturation vapor pressure over liquid
47 ! function result Real(krealfp) saturation vapor pressure in Pascals
48 ! t Real(krealfp) temperature in Kelvin
49 !
50 ! gpvsi Compute saturation vapor pressure over ice table
51 !
52 ! fpvsi Elementally compute saturation vapor pressure over ice
53 ! function result Real(krealfp) saturation vapor pressure in Pascals
54 ! t Real(krealfp) temperature in Kelvin
55 !
56 ! fpvsiq Elementally compute saturation vapor pressure over ice
57 ! function result Real(krealfp) saturation vapor pressure in Pascals
58 ! t Real(krealfp) temperature in Kelvin
59 !
60 ! fpvsix Elementally compute saturation vapor pressure over ice
61 ! function result Real(krealfp) saturation vapor pressure in Pascals
62 ! t Real(krealfp) temperature in Kelvin
63 !
64 ! gpvs Compute saturation vapor pressure table
65 !
66 ! fpvs Elementally compute saturation vapor pressure
67 ! function result Real(krealfp) saturation vapor pressure in Pascals
68 ! t Real(krealfp) temperature in Kelvin
69 !
70 ! fpvsq Elementally compute saturation vapor pressure
71 ! function result Real(krealfp) saturation vapor pressure in Pascals
72 ! t Real(krealfp) temperature in Kelvin
73 !
74 ! fpvsx Elementally compute saturation vapor pressure
75 ! function result Real(krealfp) saturation vapor pressure in Pascals
76 ! t Real(krealfp) temperature in Kelvin
77 !
78 ! gtdpl Compute dewpoint temperature over liquid table
79 !
80 ! ftdpl Elementally compute dewpoint temperature over liquid
81 ! function result Real(krealfp) dewpoint temperature in Kelvin
82 ! pv Real(krealfp) vapor pressure in Pascals
83 !
84 ! ftdplq Elementally compute dewpoint temperature over liquid
85 ! function result Real(krealfp) dewpoint temperature in Kelvin
86 ! pv Real(krealfp) vapor pressure in Pascals
87 !
88 ! ftdplx Elementally compute dewpoint temperature over liquid
89 ! function result Real(krealfp) dewpoint temperature in Kelvin
90 ! pv Real(krealfp) vapor pressure in Pascals
91 !
92 ! ftdplxg Elementally compute dewpoint temperature over liquid
93 ! function result Real(krealfp) dewpoint temperature in Kelvin
94 ! t Real(krealfp) guess dewpoint temperature in Kelvin
95 ! pv Real(krealfp) vapor pressure in Pascals
96 !
97 ! gtdpi Compute dewpoint temperature table over ice
98 !
99 ! ftdpi Elementally compute dewpoint temperature over ice
100 ! function result Real(krealfp) dewpoint temperature in Kelvin
101 ! pv Real(krealfp) vapor pressure in Pascals
102 !
103 ! ftdpiq Elementally compute dewpoint temperature over ice
104 ! function result Real(krealfp) dewpoint temperature in Kelvin
105 ! pv Real(krealfp) vapor pressure in Pascals
106 !
107 ! ftdpix Elementally compute dewpoint temperature over ice
108 ! function result Real(krealfp) dewpoint temperature in Kelvin
109 ! pv Real(krealfp) vapor pressure in Pascals
110 !
111 ! ftdpixg Elementally compute dewpoint temperature over ice
112 ! function result Real(krealfp) dewpoint temperature in Kelvin
113 ! t Real(krealfp) guess dewpoint temperature in Kelvin
114 ! pv Real(krealfp) vapor pressure in Pascals
115 !
116 ! gtdp Compute dewpoint temperature table
117 !
118 ! ftdp Elementally compute dewpoint temperature
119 ! function result Real(krealfp) dewpoint temperature in Kelvin
120 ! pv Real(krealfp) vapor pressure in Pascals
121 !
122 ! ftdpq Elementally compute dewpoint temperature
123 ! function result Real(krealfp) dewpoint temperature in Kelvin
124 ! pv Real(krealfp) vapor pressure in Pascals
125 !
126 ! ftdpx Elementally compute dewpoint temperature
127 ! function result Real(krealfp) dewpoint temperature in Kelvin
128 ! pv Real(krealfp) vapor pressure in Pascals
129 !
130 ! ftdpxg Elementally compute dewpoint temperature
131 ! function result Real(krealfp) dewpoint temperature in Kelvin
132 ! t Real(krealfp) guess dewpoint temperature in Kelvin
133 ! pv Real(krealfp) vapor pressure in Pascals
134 !
135 ! gthe Compute equivalent potential temperature table
136 !
137 ! fthe Elementally compute equivalent potential temperature
138 ! function result Real(krealfp) equivalent potential temperature in Kelvin
139 ! t Real(krealfp) LCL temperature in Kelvin
140 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
141 !
142 ! ftheq Elementally compute equivalent potential temperature
143 ! function result Real(krealfp) equivalent potential temperature in Kelvin
144 ! t Real(krealfp) LCL temperature in Kelvin
145 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
146 !
147 ! fthex Elementally compute equivalent potential temperature
148 ! function result Real(krealfp) equivalent potential temperature in Kelvin
149 ! t Real(krealfp) LCL temperature in Kelvin
150 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
151 !
152 ! gtma Compute moist adiabat tables
153 !
154 ! stma Elementally compute moist adiabat temperature and moisture
155 ! the Real(krealfp) equivalent potential temperature in Kelvin
156 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
157 ! tma Real(krealfp) parcel temperature in Kelvin
158 ! qma Real(krealfp) parcel specific humidity in kg/kg
159 !
160 ! stmaq Elementally compute moist adiabat temperature and moisture
161 ! the Real(krealfp) equivalent potential temperature in Kelvin
162 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
163 ! tma Real(krealfp) parcel temperature in Kelvin
164 ! qma Real(krealfp) parcel specific humidity in kg/kg
165 !
166 ! stmax Elementally compute moist adiabat temperature and moisture
167 ! the Real(krealfp) equivalent potential temperature in Kelvin
168 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
169 ! tma Real(krealfp) parcel temperature in Kelvin
170 ! qma Real(krealfp) parcel specific humidity in kg/kg
171 !
172 ! stmaxg Elementally compute moist adiabat temperature and moisture
173 ! tg Real(krealfp) guess parcel temperature in Kelvin
174 ! the Real(krealfp) equivalent potential temperature in Kelvin
175 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
176 ! tma Real(krealfp) parcel temperature in Kelvin
177 ! qma Real(krealfp) parcel specific humidity in kg/kg
178 !
179 ! gpkap Compute pressure to the kappa table
180 !
181 ! fpkap Elementally raise pressure to the kappa power.
182 ! function result Real(krealfp) p over 1e5 Pa to the kappa power
183 ! p Real(krealfp) pressure in Pascals
184 !
185 ! fpkapq Elementally raise pressure to the kappa power.
186 ! function result Real(krealfp) p over 1e5 Pa to the kappa power
187 ! p Real(krealfp) pressure in Pascals
188 !
189 ! fpkapo Elementally raise pressure to the kappa power.
190 ! function result Real(krealfp) p over 1e5 Pa to the kappa power
191 ! p Real(krealfp) surface pressure in Pascals
192 !
193 ! fpkapx Elementally raise pressure to the kappa power.
194 ! function result Real(krealfp) p over 1e5 Pa to the kappa power
195 ! p Real(krealfp) pressure in Pascals
196 !
197 ! grkap Compute pressure to the 1/kappa table
198 !
199 ! frkap Elementally raise pressure to the 1/kappa power.
200 ! function result Real(krealfp) pressure in Pascals
201 ! pkap Real(krealfp) p over 1e5 Pa to the 1/kappa power
202 !
203 ! frkapq Elementally raise pressure to the kappa power.
204 ! function result Real(krealfp) pressure in Pascals
205 ! pkap Real(krealfp) p over 1e5 Pa to the kappa power
206 !
207 ! frkapx Elementally raise pressure to the kappa power.
208 ! function result Real(krealfp) pressure in Pascals
209 ! pkap Real(krealfp) p over 1e5 Pa to the kappa power
210 !
211 ! gtlcl Compute LCL temperature table
212 !
213 ! ftlcl Elementally compute LCL temperature.
214 ! function result Real(krealfp) temperature at the LCL in Kelvin
215 ! t Real(krealfp) temperature in Kelvin
216 ! tdpd Real(krealfp) dewpoint depression in Kelvin
217 !
218 ! ftlclq Elementally compute LCL temperature.
219 ! function result Real(krealfp) temperature at the LCL in Kelvin
220 ! t Real(krealfp) temperature in Kelvin
221 ! tdpd Real(krealfp) dewpoint depression in Kelvin
222 !
223 ! ftlclo Elementally compute LCL temperature.
224 ! function result Real(krealfp) temperature at the LCL in Kelvin
225 ! t Real(krealfp) temperature in Kelvin
226 ! tdpd Real(krealfp) dewpoint depression in Kelvin
227 !
228 ! ftlclx Elementally compute LCL temperature.
229 ! function result Real(krealfp) temperature at the LCL in Kelvin
230 ! t Real(krealfp) temperature in Kelvin
231 ! tdpd Real(krealfp) dewpoint depression in Kelvin
232 !
233 ! gfuncphys Compute all physics function tables
234 !
235 ! Attributes:
236 ! Language: Fortran 90
237 !
238 !$$$
239  use machine,only:kind_phys
240  use physcons
241  implicit none
242  private
243 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244 ! Public Variables
245 ! integer,public,parameter:: krealfp=selected_real_kind(15,45)
246  integer,public,parameter:: krealfp=kind_phys
247 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248 ! Private Variables
249  real(krealfp),parameter:: psatb=con_psat*1.e-5
250  integer,parameter:: nxpvsl=7501
251  real(krealfp) c1xpvsl,c2xpvsl,tbpvsl(nxpvsl)
252  integer,parameter:: nxpvsi=7501
253  real(krealfp) c1xpvsi,c2xpvsi,tbpvsi(nxpvsi)
254  integer,parameter:: nxpvs=7501
255  real(krealfp) c1xpvs,c2xpvs,tbpvs(nxpvs)
256  integer,parameter:: nxtdpl=5001
257  real(krealfp) c1xtdpl,c2xtdpl,tbtdpl(nxtdpl)
258  integer,parameter:: nxtdpi=5001
259  real(krealfp) c1xtdpi,c2xtdpi,tbtdpi(nxtdpi)
260  integer,parameter:: nxtdp=5001
261  real(krealfp) c1xtdp,c2xtdp,tbtdp(nxtdp)
262  integer,parameter:: nxthe=241,nythe=151
264  integer,parameter:: nxma=151,nyma=121
266 ! integer,parameter:: nxpkap=5501
267  integer,parameter:: nxpkap=11001
268  real(krealfp) c1xpkap,c2xpkap,tbpkap(nxpkap)
269  integer,parameter:: nxrkap=5501
270  real(krealfp) c1xrkap,c2xrkap,tbrkap(nxrkap)
271  integer,parameter:: nxtlcl=151,nytlcl=61
273 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274 ! Public Subprograms
275  public gpvsl,fpvsl,fpvslq,fpvslx
276  public gpvsi,fpvsi,fpvsiq,fpvsix
277  public gpvs,fpvs,fpvsq,fpvsx
280  public gtdp,ftdp,ftdpq,ftdpx,ftdpxg
281  public gthe,fthe,ftheq,fthex
282  public gtma,stma,stmaq,stmax,stmaxg
284  public grkap,frkap,frkapq,frkapx
286  public gfuncphys
287 contains
288 !-------------------------------------------------------------------------------
289  subroutine gpvsl
290 !$$$ Subprogram Documentation Block
291 !
292 ! Subprogram: gpvsl Compute saturation vapor pressure table over liquid
293 ! Author: N Phillips W/NMC2X2 Date: 30 dec 82
294 !
295 ! Abstract: Computes saturation vapor pressure table as a function of
296 ! temperature for the table lookup function fpvsl.
297 ! Exact saturation vapor pressures are calculated in subprogram fpvslx.
298 ! The current implementation computes a table with a length
299 ! of 7501 for temperatures ranging from 180. to 330. Kelvin.
300 !
301 ! Program History Log:
302 ! 91-05-07 Iredell
303 ! 94-12-30 Iredell expand table
304 ! 1999-03-01 Iredell f90 module
305 !
306 ! Usage: call gpvsl
307 !
308 ! Subprograms called:
309 ! (fpvslx) inlinable function to compute saturation vapor pressure
310 !
311 ! Attributes:
312 ! Language: Fortran 90.
313 !
314 !$$$
315  implicit none
316  integer jx
317  real(krealfp) xmin,xmax,xinc,x,t
318 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
319  xmin=180.0_krealfp
320  xmax=330.0_krealfp
321  xinc=(xmax-xmin)/(nxpvsl-1)
322 ! c1xpvsl=1.-xmin/xinc
323  c2xpvsl=1./xinc
324  c1xpvsl=1.-xmin*c2xpvsl
325  do jx=1,nxpvsl
326  x=xmin+(jx-1)*xinc
327  t=x
328  tbpvsl(jx)=fpvslx(t)
329  enddo
330 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
331  end subroutine
332 !-------------------------------------------------------------------------------
333  elemental function fpvsl(t)
334 !$$$ Subprogram Documentation Block
335 !
336 ! Subprogram: fpvsl Compute saturation vapor pressure over liquid
337 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
338 !
339 ! Abstract: Compute saturation vapor pressure from the temperature.
340 ! A linear interpolation is done between values in a lookup table
341 ! computed in gpvsl. See documentation for fpvslx for details.
342 ! Input values outside table range are reset to table extrema.
343 ! The interpolation accuracy is almost 6 decimal places.
344 ! On the Cray, fpvsl is about 4 times faster than exact calculation.
345 ! This function should be expanded inline in the calling routine.
346 !
347 ! Program History Log:
348 ! 91-05-07 Iredell made into inlinable function
349 ! 94-12-30 Iredell expand table
350 ! 1999-03-01 Iredell f90 module
351 !
352 ! Usage: pvsl=fpvsl(t)
353 !
354 ! Input argument list:
355 ! t Real(krealfp) temperature in Kelvin
356 !
357 ! Output argument list:
358 ! fpvsl Real(krealfp) saturation vapor pressure in Pascals
359 !
360 ! Attributes:
361 ! Language: Fortran 90.
362 !
363 !$$$
364  implicit none
365  real(krealfp) fpvsl
366  real(krealfp),intent(in):: t
367  integer jx
368  real(krealfp) xj
369 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
370  xj=min(max(c1xpvsl+c2xpvsl*t,1._krealfp),real(nxpvsl,krealfp))
371  jx=min(xj,nxpvsl-1._krealfp)
372  fpvsl=tbpvsl(jx)+(xj-jx)*(tbpvsl(jx+1)-tbpvsl(jx))
373 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
374  end function
375 !-------------------------------------------------------------------------------
376  elemental function fpvslq(t)
377 !$$$ Subprogram Documentation Block
378 !
379 ! Subprogram: fpvslq Compute saturation vapor pressure over liquid
380 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
381 !
382 ! Abstract: Compute saturation vapor pressure from the temperature.
383 ! A quadratic interpolation is done between values in a lookup table
384 ! computed in gpvsl. See documentation for fpvslx for details.
385 ! Input values outside table range are reset to table extrema.
386 ! The interpolation accuracy is almost 9 decimal places.
387 ! On the Cray, fpvslq is about 3 times faster than exact calculation.
388 ! This function should be expanded inline in the calling routine.
389 !
390 ! Program History Log:
391 ! 91-05-07 Iredell made into inlinable function
392 ! 94-12-30 Iredell quadratic interpolation
393 ! 1999-03-01 Iredell f90 module
394 !
395 ! Usage: pvsl=fpvslq(t)
396 !
397 ! Input argument list:
398 ! t Real(krealfp) temperature in Kelvin
399 !
400 ! Output argument list:
401 ! fpvslq Real(krealfp) saturation vapor pressure in Pascals
402 !
403 ! Attributes:
404 ! Language: Fortran 90.
405 !
406 !$$$
407  implicit none
408  real(krealfp) fpvslq
409  real(krealfp),intent(in):: t
410  integer jx
411  real(krealfp) xj,dxj,fj1,fj2,fj3
412 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
413  xj=min(max(c1xpvsl+c2xpvsl*t,1._krealfp),real(nxpvsl,krealfp))
414  jx=min(max(nint(xj),2),nxpvsl-1)
415  dxj=xj-jx
416  fj1=tbpvsl(jx-1)
417  fj2=tbpvsl(jx)
418  fj3=tbpvsl(jx+1)
419  fpvslq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
420 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
421  end function
422 !-------------------------------------------------------------------------------
423  elemental function fpvslx(t)
424 !$$$ Subprogram Documentation Block
425 !
426 ! Subprogram: fpvslx Compute saturation vapor pressure over liquid
427 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
428 !
429 ! Abstract: Exactly compute saturation vapor pressure from temperature.
430 ! The water model assumes a perfect gas, constant specific heats
431 ! for gas and liquid, and neglects the volume of the liquid.
432 ! The model does account for the variation of the latent heat
433 ! of condensation with temperature. The ice option is not included.
434 ! The Clausius-Clapeyron equation is integrated from the triple point
435 ! to get the formula
436 ! pvsl=con_psat*(tr**xa)*exp(xb*(1.-tr))
437 ! where tr is ttp/t and other values are physical constants.
438 ! This function should be expanded inline in the calling routine.
439 !
440 ! Program History Log:
441 ! 91-05-07 Iredell made into inlinable function
442 ! 94-12-30 Iredell exact computation
443 ! 1999-03-01 Iredell f90 module
444 !
445 ! Usage: pvsl=fpvslx(t)
446 !
447 ! Input argument list:
448 ! t Real(krealfp) temperature in Kelvin
449 !
450 ! Output argument list:
451 ! fpvslx Real(krealfp) saturation vapor pressure in Pascals
452 !
453 ! Attributes:
454 ! Language: Fortran 90.
455 !
456 !$$$
457  implicit none
458  real(krealfp) fpvslx
459  real(krealfp),intent(in):: t
460  real(krealfp),parameter:: dldt=con_cvap-con_cliq
461  real(krealfp),parameter:: heat=con_hvap
462  real(krealfp),parameter:: xpona=-dldt/con_rv
463  real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
464  real(krealfp) tr
465 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
466  tr=con_ttp/t
467  fpvslx=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
468 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
469  end function
470 !-------------------------------------------------------------------------------
471  subroutine gpvsi
472 !$$$ Subprogram Documentation Block
473 !
474 ! Subprogram: gpvsi Compute saturation vapor pressure table over ice
475 ! Author: N Phillips W/NMC2X2 Date: 30 dec 82
476 !
477 ! Abstract: Computes saturation vapor pressure table as a function of
478 ! temperature for the table lookup function fpvsi.
479 ! Exact saturation vapor pressures are calculated in subprogram fpvsix.
480 ! The current implementation computes a table with a length
481 ! of 7501 for temperatures ranging from 180. to 330. Kelvin.
482 !
483 ! Program History Log:
484 ! 91-05-07 Iredell
485 ! 94-12-30 Iredell expand table
486 ! 1999-03-01 Iredell f90 module
487 ! 2001-02-26 Iredell ice phase
488 !
489 ! Usage: call gpvsi
490 !
491 ! Subprograms called:
492 ! (fpvsix) inlinable function to compute saturation vapor pressure
493 !
494 ! Attributes:
495 ! Language: Fortran 90.
496 !
497 !$$$
498  implicit none
499  integer jx
500  real(krealfp) xmin,xmax,xinc,x,t
501 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
502  xmin=180.0_krealfp
503  xmax=330.0_krealfp
504  xinc=(xmax-xmin)/(nxpvsi-1)
505 ! c1xpvsi=1.-xmin/xinc
506  c2xpvsi=1./xinc
507  c1xpvsi=1.-xmin*c2xpvsi
508  do jx=1,nxpvsi
509  x=xmin+(jx-1)*xinc
510  t=x
511  tbpvsi(jx)=fpvsix(t)
512  enddo
513 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
514  end subroutine
515 !-------------------------------------------------------------------------------
516  elemental function fpvsi(t)
517 !$$$ Subprogram Documentation Block
518 !
519 ! Subprogram: fpvsi Compute saturation vapor pressure over ice
520 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
521 !
522 ! Abstract: Compute saturation vapor pressure from the temperature.
523 ! A linear interpolation is done between values in a lookup table
524 ! computed in gpvsi. See documentation for fpvsix for details.
525 ! Input values outside table range are reset to table extrema.
526 ! The interpolation accuracy is almost 6 decimal places.
527 ! On the Cray, fpvsi is about 4 times faster than exact calculation.
528 ! This function should be expanded inline in the calling routine.
529 !
530 ! Program History Log:
531 ! 91-05-07 Iredell made into inlinable function
532 ! 94-12-30 Iredell expand table
533 ! 1999-03-01 Iredell f90 module
534 ! 2001-02-26 Iredell ice phase
535 !
536 ! Usage: pvsi=fpvsi(t)
537 !
538 ! Input argument list:
539 ! t Real(krealfp) temperature in Kelvin
540 !
541 ! Output argument list:
542 ! fpvsi Real(krealfp) saturation vapor pressure in Pascals
543 !
544 ! Attributes:
545 ! Language: Fortran 90.
546 !
547 !$$$
548  implicit none
549  real(krealfp) fpvsi
550  real(krealfp),intent(in):: t
551  integer jx
552  real(krealfp) xj
553 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
554  xj=min(max(c1xpvsi+c2xpvsi*t,1._krealfp),real(nxpvsi,krealfp))
555  jx=min(xj,nxpvsi-1._krealfp)
556  fpvsi=tbpvsi(jx)+(xj-jx)*(tbpvsi(jx+1)-tbpvsi(jx))
557 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
558  end function
559 !-------------------------------------------------------------------------------
560  elemental function fpvsiq(t)
561 !$$$ Subprogram Documentation Block
562 !
563 ! Subprogram: fpvsiq Compute saturation vapor pressure over ice
564 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
565 !
566 ! Abstract: Compute saturation vapor pressure from the temperature.
567 ! A quadratic interpolation is done between values in a lookup table
568 ! computed in gpvsi. See documentation for fpvsix for details.
569 ! Input values outside table range are reset to table extrema.
570 ! The interpolation accuracy is almost 9 decimal places.
571 ! On the Cray, fpvsiq is about 3 times faster than exact calculation.
572 ! This function should be expanded inline in the calling routine.
573 !
574 ! Program History Log:
575 ! 91-05-07 Iredell made into inlinable function
576 ! 94-12-30 Iredell quadratic interpolation
577 ! 1999-03-01 Iredell f90 module
578 ! 2001-02-26 Iredell ice phase
579 !
580 ! Usage: pvsi=fpvsiq(t)
581 !
582 ! Input argument list:
583 ! t Real(krealfp) temperature in Kelvin
584 !
585 ! Output argument list:
586 ! fpvsiq Real(krealfp) saturation vapor pressure in Pascals
587 !
588 ! Attributes:
589 ! Language: Fortran 90.
590 !
591 !$$$
592  implicit none
593  real(krealfp) fpvsiq
594  real(krealfp),intent(in):: t
595  integer jx
596  real(krealfp) xj,dxj,fj1,fj2,fj3
597 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
598  xj=min(max(c1xpvsi+c2xpvsi*t,1._krealfp),real(nxpvsi,krealfp))
599  jx=min(max(nint(xj),2),nxpvsi-1)
600  dxj=xj-jx
601  fj1=tbpvsi(jx-1)
602  fj2=tbpvsi(jx)
603  fj3=tbpvsi(jx+1)
604  fpvsiq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
605 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
606  end function
607 !-------------------------------------------------------------------------------
608  elemental function fpvsix(t)
609 !$$$ Subprogram Documentation Block
610 !
611 ! Subprogram: fpvsix Compute saturation vapor pressure over ice
612 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
613 !
614 ! Abstract: Exactly compute saturation vapor pressure from temperature.
615 ! The water model assumes a perfect gas, constant specific heats
616 ! for gas and ice, and neglects the volume of the ice.
617 ! The model does account for the variation of the latent heat
618 ! of condensation with temperature. The liquid option is not included.
619 ! The Clausius-Clapeyron equation is integrated from the triple point
620 ! to get the formula
621 ! pvsi=con_psat*(tr**xa)*exp(xb*(1.-tr))
622 ! where tr is ttp/t and other values are physical constants.
623 ! This function should be expanded inline in the calling routine.
624 !
625 ! Program History Log:
626 ! 91-05-07 Iredell made into inlinable function
627 ! 94-12-30 Iredell exact computation
628 ! 1999-03-01 Iredell f90 module
629 ! 2001-02-26 Iredell ice phase
630 !
631 ! Usage: pvsi=fpvsix(t)
632 !
633 ! Input argument list:
634 ! t Real(krealfp) temperature in Kelvin
635 !
636 ! Output argument list:
637 ! fpvsix Real(krealfp) saturation vapor pressure in Pascals
638 !
639 ! Attributes:
640 ! Language: Fortran 90.
641 !
642 !$$$
643  implicit none
644  real(krealfp) fpvsix
645  real(krealfp),intent(in):: t
646  real(krealfp),parameter:: dldt=con_cvap-con_csol
647  real(krealfp),parameter:: heat=con_hvap+con_hfus
648  real(krealfp),parameter:: xpona=-dldt/con_rv
649  real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
650  real(krealfp) tr
651 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
652  tr=con_ttp/t
653  fpvsix=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
654 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
655  end function
656 !-------------------------------------------------------------------------------
657  subroutine gpvs
658 !$$$ Subprogram Documentation Block
659 !
660 ! Subprogram: gpvs Compute saturation vapor pressure table
661 ! Author: N Phillips W/NMC2X2 Date: 30 dec 82
662 !
663 ! Abstract: Computes saturation vapor pressure table as a function of
664 ! temperature for the table lookup function fpvs.
665 ! Exact saturation vapor pressures are calculated in subprogram fpvsx.
666 ! The current implementation computes a table with a length
667 ! of 7501 for temperatures ranging from 180. to 330. Kelvin.
668 !
669 ! Program History Log:
670 ! 91-05-07 Iredell
671 ! 94-12-30 Iredell expand table
672 ! 1999-03-01 Iredell f90 module
673 ! 2001-02-26 Iredell ice phase
674 !
675 ! Usage: call gpvs
676 !
677 ! Subprograms called:
678 ! (fpvsx) inlinable function to compute saturation vapor pressure
679 !
680 ! Attributes:
681 ! Language: Fortran 90.
682 !
683 !$$$
684  implicit none
685  integer jx
686  real(krealfp) xmin,xmax,xinc,x,t
687 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
688  xmin=180.0_krealfp
689  xmax=330.0_krealfp
690  xinc=(xmax-xmin)/(nxpvs-1)
691 ! c1xpvs=1.-xmin/xinc
692  c2xpvs=1./xinc
693  c1xpvs=1.-xmin*c2xpvs
694  do jx=1,nxpvs
695  x=xmin+(jx-1)*xinc
696  t=x
697  tbpvs(jx)=fpvsx(t)
698  enddo
699 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
700  end subroutine
701 !-------------------------------------------------------------------------------
702  elemental function fpvs(t)
703 !$$$ Subprogram Documentation Block
704 !
705 ! Subprogram: fpvs Compute saturation vapor pressure
706 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
707 !
708 ! Abstract: Compute saturation vapor pressure from the temperature.
709 ! A linear interpolation is done between values in a lookup table
710 ! computed in gpvs. See documentation for fpvsx for details.
711 ! Input values outside table range are reset to table extrema.
712 ! The interpolation accuracy is almost 6 decimal places.
713 ! On the Cray, fpvs is about 4 times faster than exact calculation.
714 ! This function should be expanded inline in the calling routine.
715 !
716 ! Program History Log:
717 ! 91-05-07 Iredell made into inlinable function
718 ! 94-12-30 Iredell expand table
719 ! 1999-03-01 Iredell f90 module
720 ! 2001-02-26 Iredell ice phase
721 !
722 ! Usage: pvs=fpvs(t)
723 !
724 ! Input argument list:
725 ! t Real(krealfp) temperature in Kelvin
726 !
727 ! Output argument list:
728 ! fpvs Real(krealfp) saturation vapor pressure in Pascals
729 !
730 ! Attributes:
731 ! Language: Fortran 90.
732 !
733 !$$$
734  implicit none
735  real(krealfp) fpvs
736  real(krealfp),intent(in):: t
737  integer jx
738  real(krealfp) xj
739 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
740  xj=min(max(c1xpvs+c2xpvs*t,1._krealfp),real(nxpvs,krealfp))
741  jx=min(xj,nxpvs-1._krealfp)
742  fpvs=tbpvs(jx)+(xj-jx)*(tbpvs(jx+1)-tbpvs(jx))
743 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
744  end function
745 !-------------------------------------------------------------------------------
746  elemental function fpvsq(t)
747 !$$$ Subprogram Documentation Block
748 !
749 ! Subprogram: fpvsq Compute saturation vapor pressure
750 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
751 !
752 ! Abstract: Compute saturation vapor pressure from the temperature.
753 ! A quadratic interpolation is done between values in a lookup table
754 ! computed in gpvs. See documentation for fpvsx for details.
755 ! Input values outside table range are reset to table extrema.
756 ! The interpolation accuracy is almost 9 decimal places.
757 ! On the Cray, fpvsq is about 3 times faster than exact calculation.
758 ! This function should be expanded inline in the calling routine.
759 !
760 ! Program History Log:
761 ! 91-05-07 Iredell made into inlinable function
762 ! 94-12-30 Iredell quadratic interpolation
763 ! 1999-03-01 Iredell f90 module
764 ! 2001-02-26 Iredell ice phase
765 !
766 ! Usage: pvs=fpvsq(t)
767 !
768 ! Input argument list:
769 ! t Real(krealfp) temperature in Kelvin
770 !
771 ! Output argument list:
772 ! fpvsq Real(krealfp) saturation vapor pressure in Pascals
773 !
774 ! Attributes:
775 ! Language: Fortran 90.
776 !
777 !$$$
778  implicit none
779  real(krealfp) fpvsq
780  real(krealfp),intent(in):: t
781  integer jx
782  real(krealfp) xj,dxj,fj1,fj2,fj3
783 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
784  xj=min(max(c1xpvs+c2xpvs*t,1._krealfp),real(nxpvs,krealfp))
785  jx=min(max(nint(xj),2),nxpvs-1)
786  dxj=xj-jx
787  fj1=tbpvs(jx-1)
788  fj2=tbpvs(jx)
789  fj3=tbpvs(jx+1)
790  fpvsq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
791 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
792  end function
793 !-------------------------------------------------------------------------------
794  elemental function fpvsx(t)
795 !$$$ Subprogram Documentation Block
796 !
797 ! Subprogram: fpvsx Compute saturation vapor pressure
798 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
799 !
800 ! Abstract: Exactly compute saturation vapor pressure from temperature.
801 ! The saturation vapor pressure over either liquid and ice is computed
802 ! over liquid for temperatures above the triple point,
803 ! over ice for temperatures 20 degress below the triple point,
804 ! and a linear combination of the two for temperatures in between.
805 ! The water model assumes a perfect gas, constant specific heats
806 ! for gas, liquid and ice, and neglects the volume of the condensate.
807 ! The model does account for the variation of the latent heat
808 ! of condensation and sublimation with temperature.
809 ! The Clausius-Clapeyron equation is integrated from the triple point
810 ! to get the formula
811 ! pvsl=con_psat*(tr**xa)*exp(xb*(1.-tr))
812 ! where tr is ttp/t and other values are physical constants.
813 ! The reference for this computation is Emanuel(1994), pages 116-117.
814 ! This function should be expanded inline in the calling routine.
815 !
816 ! Program History Log:
817 ! 91-05-07 Iredell made into inlinable function
818 ! 94-12-30 Iredell exact computation
819 ! 1999-03-01 Iredell f90 module
820 ! 2001-02-26 Iredell ice phase
821 !
822 ! Usage: pvs=fpvsx(t)
823 !
824 ! Input argument list:
825 ! t Real(krealfp) temperature in Kelvin
826 !
827 ! Output argument list:
828 ! fpvsx Real(krealfp) saturation vapor pressure in Pascals
829 !
830 ! Attributes:
831 ! Language: Fortran 90.
832 !
833 !$$$
834  implicit none
835  real(krealfp) fpvsx
836  real(krealfp),intent(in):: t
837  real(krealfp),parameter:: tliq=con_ttp
838  real(krealfp),parameter:: tice=con_ttp-20.0
839  real(krealfp),parameter:: dldtl=con_cvap-con_cliq
840  real(krealfp),parameter:: heatl=con_hvap
841  real(krealfp),parameter:: xponal=-dldtl/con_rv
842  real(krealfp),parameter:: xponbl=-dldtl/con_rv+heatl/(con_rv*con_ttp)
843  real(krealfp),parameter:: dldti=con_cvap-con_csol
844  real(krealfp),parameter:: heati=con_hvap+con_hfus
845  real(krealfp),parameter:: xponai=-dldti/con_rv
846  real(krealfp),parameter:: xponbi=-dldti/con_rv+heati/(con_rv*con_ttp)
847  real(krealfp) tr,w,pvl,pvi
848 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
849  tr=con_ttp/t
850  if(t.ge.tliq) then
851  fpvsx=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
852  elseif(t.lt.tice) then
853  fpvsx=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
854  else
855  w=(t-tice)/(tliq-tice)
856  pvl=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
857  pvi=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
858  fpvsx=w*pvl+(1.-w)*pvi
859  endif
860 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
861  end function
862 !-------------------------------------------------------------------------------
863  subroutine gtdpl
864 !$$$ Subprogram Documentation Block
865 !
866 ! Subprogram: gtdpl Compute dewpoint temperature over liquid table
867 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
868 !
869 ! Abstract: Compute dewpoint temperature table as a function of
870 ! vapor pressure for inlinable function ftdpl.
871 ! Exact dewpoint temperatures are calculated in subprogram ftdplxg.
872 ! The current implementation computes a table with a length
873 ! of 5001 for vapor pressures ranging from 1 to 10001 Pascals
874 ! giving a dewpoint temperature range of 208 to 319 Kelvin.
875 !
876 ! Program History Log:
877 ! 91-05-07 Iredell
878 ! 94-12-30 Iredell expand table
879 ! 1999-03-01 Iredell f90 module
880 !
881 ! Usage: call gtdpl
882 !
883 ! Subprograms called:
884 ! (ftdplxg) inlinable function to compute dewpoint temperature over liquid
885 !
886 ! Attributes:
887 ! Language: Fortran 90.
888 !
889 !$$$
890  implicit none
891  integer jx
892  real(krealfp) xmin,xmax,xinc,t,x,pv
893 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
894  xmin=1
895  xmax=10001
896  xinc=(xmax-xmin)/(nxtdpl-1)
897  c1xtdpl=1.-xmin/xinc
898  c2xtdpl=1./xinc
899  t=208.0
900  do jx=1,nxtdpl
901  x=xmin+(jx-1)*xinc
902  pv=x
903  t=ftdplxg(t,pv)
904  tbtdpl(jx)=t
905  enddo
906 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
907  end subroutine
908 !-------------------------------------------------------------------------------
909  elemental function ftdpl(pv)
910 !$$$ Subprogram Documentation Block
911 !
912 ! Subprogram: ftdpl Compute dewpoint temperature over liquid
913 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
914 !
915 ! Abstract: Compute dewpoint temperature from vapor pressure.
916 ! A linear interpolation is done between values in a lookup table
917 ! computed in gtdpl. See documentation for ftdplxg for details.
918 ! Input values outside table range are reset to table extrema.
919 ! The interpolation accuracy is better than 0.0005 Kelvin
920 ! for dewpoint temperatures greater than 250 Kelvin,
921 ! but decreases to 0.02 Kelvin for a dewpoint around 230 Kelvin.
922 ! On the Cray, ftdpl is about 75 times faster than exact calculation.
923 ! This function should be expanded inline in the calling routine.
924 !
925 ! Program History Log:
926 ! 91-05-07 Iredell made into inlinable function
927 ! 94-12-30 Iredell expand table
928 ! 1999-03-01 Iredell f90 module
929 !
930 ! Usage: tdpl=ftdpl(pv)
931 !
932 ! Input argument list:
933 ! pv Real(krealfp) vapor pressure in Pascals
934 !
935 ! Output argument list:
936 ! ftdpl Real(krealfp) dewpoint temperature in Kelvin
937 !
938 ! Attributes:
939 ! Language: Fortran 90.
940 !
941 !$$$
942  implicit none
943  real(krealfp) ftdpl
944  real(krealfp),intent(in):: pv
945  integer jx
946  real(krealfp) xj
947 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
948  xj=min(max(c1xtdpl+c2xtdpl*pv,1._krealfp),real(nxtdpl,krealfp))
949  jx=min(xj,nxtdpl-1._krealfp)
950  ftdpl=tbtdpl(jx)+(xj-jx)*(tbtdpl(jx+1)-tbtdpl(jx))
951 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
952  end function
953 !-------------------------------------------------------------------------------
954  elemental function ftdplq(pv)
955 !$$$ Subprogram Documentation Block
956 !
957 ! Subprogram: ftdplq Compute dewpoint temperature over liquid
958 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
959 !
960 ! Abstract: Compute dewpoint temperature from vapor pressure.
961 ! A quadratic interpolation is done between values in a lookup table
962 ! computed in gtdpl. see documentation for ftdplxg for details.
963 ! Input values outside table range are reset to table extrema.
964 ! the interpolation accuracy is better than 0.00001 Kelvin
965 ! for dewpoint temperatures greater than 250 Kelvin,
966 ! but decreases to 0.002 Kelvin for a dewpoint around 230 Kelvin.
967 ! On the Cray, ftdplq is about 60 times faster than exact calculation.
968 ! This function should be expanded inline in the calling routine.
969 !
970 ! Program History Log:
971 ! 91-05-07 Iredell made into inlinable function
972 ! 94-12-30 Iredell quadratic interpolation
973 ! 1999-03-01 Iredell f90 module
974 !
975 ! Usage: tdpl=ftdplq(pv)
976 !
977 ! Input argument list:
978 ! pv Real(krealfp) vapor pressure in Pascals
979 !
980 ! Output argument list:
981 ! ftdplq Real(krealfp) dewpoint temperature in Kelvin
982 !
983 ! Attributes:
984 ! Language: Fortran 90.
985 !
986 !$$$
987  implicit none
988  real(krealfp) ftdplq
989  real(krealfp),intent(in):: pv
990  integer jx
991  real(krealfp) xj,dxj,fj1,fj2,fj3
992 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
993  xj=min(max(c1xtdpl+c2xtdpl*pv,1._krealfp),real(nxtdpl,krealfp))
994  jx=min(max(nint(xj),2),nxtdpl-1)
995  dxj=xj-jx
996  fj1=tbtdpl(jx-1)
997  fj2=tbtdpl(jx)
998  fj3=tbtdpl(jx+1)
999  ftdplq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
1000 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1001  end function
1002 !-------------------------------------------------------------------------------
1003  elemental function ftdplx(pv)
1004 !$$$ Subprogram Documentation Block
1005 !
1006 ! Subprogram: ftdplx Compute dewpoint temperature over liquid
1007 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1008 !
1009 ! Abstract: exactly compute dewpoint temperature from vapor pressure.
1010 ! An approximate dewpoint temperature for function ftdplxg
1011 ! is obtained using ftdpl so gtdpl must be already called.
1012 ! See documentation for ftdplxg for details.
1013 !
1014 ! Program History Log:
1015 ! 91-05-07 Iredell made into inlinable function
1016 ! 94-12-30 Iredell exact computation
1017 ! 1999-03-01 Iredell f90 module
1018 !
1019 ! Usage: tdpl=ftdplx(pv)
1020 !
1021 ! Input argument list:
1022 ! pv Real(krealfp) vapor pressure in Pascals
1023 !
1024 ! Output argument list:
1025 ! ftdplx Real(krealfp) dewpoint temperature in Kelvin
1026 !
1027 ! Subprograms called:
1028 ! (ftdpl) inlinable function to compute dewpoint temperature over liquid
1029 ! (ftdplxg) inlinable function to compute dewpoint temperature over liquid
1030 !
1031 ! Attributes:
1032 ! Language: Fortran 90.
1033 !
1034 !$$$
1035  implicit none
1036  real(krealfp) ftdplx
1037  real(krealfp),intent(in):: pv
1038  real(krealfp) tg
1039 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1040  tg=ftdpl(pv)
1041  ftdplx=ftdplxg(tg,pv)
1042 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1043  end function
1044 !-------------------------------------------------------------------------------
1045  elemental function ftdplxg(tg,pv)
1046 !$$$ Subprogram Documentation Block
1047 !
1048 ! Subprogram: ftdplxg Compute dewpoint temperature over liquid
1049 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1050 !
1051 ! Abstract: Exactly compute dewpoint temperature from vapor pressure.
1052 ! A guess dewpoint temperature must be provided.
1053 ! The water model assumes a perfect gas, constant specific heats
1054 ! for gas and liquid, and neglects the volume of the liquid.
1055 ! The model does account for the variation of the latent heat
1056 ! of condensation with temperature. The ice option is not included.
1057 ! The Clausius-Clapeyron equation is integrated from the triple point
1058 ! to get the formula
1059 ! pvs=con_psat*(tr**xa)*exp(xb*(1.-tr))
1060 ! where tr is ttp/t and other values are physical constants.
1061 ! The formula is inverted by iterating Newtonian approximations
1062 ! for each pvs until t is found to within 1.e-6 Kelvin.
1063 ! This function can be expanded inline in the calling routine.
1064 !
1065 ! Program History Log:
1066 ! 91-05-07 Iredell made into inlinable function
1067 ! 94-12-30 Iredell exact computation
1068 ! 1999-03-01 Iredell f90 module
1069 !
1070 ! Usage: tdpl=ftdplxg(tg,pv)
1071 !
1072 ! Input argument list:
1073 ! tg Real(krealfp) guess dewpoint temperature in Kelvin
1074 ! pv Real(krealfp) vapor pressure in Pascals
1075 !
1076 ! Output argument list:
1077 ! ftdplxg Real(krealfp) dewpoint temperature in Kelvin
1078 !
1079 ! Attributes:
1080 ! Language: Fortran 90.
1081 !
1082 !$$$
1083  implicit none
1084  real(krealfp) ftdplxg
1085  real(krealfp),intent(in):: tg,pv
1086  real(krealfp),parameter:: terrm=1.e-6
1087  real(krealfp),parameter:: dldt=con_cvap-con_cliq
1088  real(krealfp),parameter:: heat=con_hvap
1089  real(krealfp),parameter:: xpona=-dldt/con_rv
1090  real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
1091  real(krealfp) t,tr,pvt,el,dpvt,terr
1092  integer i
1093 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1094  t=tg
1095  do i=1,100
1096  tr=con_ttp/t
1097  pvt=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
1098  el=heat+dldt*(t-con_ttp)
1099  dpvt=el*pvt/(con_rv*t**2)
1100  terr=(pvt-pv)/dpvt
1101  t=t-terr
1102  if(abs(terr).le.terrm) exit
1103  enddo
1104  ftdplxg=t
1105 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1106  end function
1107 !-------------------------------------------------------------------------------
1108  subroutine gtdpi
1109 !$$$ Subprogram Documentation Block
1110 !
1111 ! Subprogram: gtdpi Compute dewpoint temperature over ice table
1112 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1113 !
1114 ! Abstract: Compute dewpoint temperature table as a function of
1115 ! vapor pressure for inlinable function ftdpi.
1116 ! Exact dewpoint temperatures are calculated in subprogram ftdpixg.
1117 ! The current implementation computes a table with a length
1118 ! of 5001 for vapor pressures ranging from 0.1 to 1000.1 Pascals
1119 ! giving a dewpoint temperature range of 197 to 279 Kelvin.
1120 !
1121 ! Program History Log:
1122 ! 91-05-07 Iredell
1123 ! 94-12-30 Iredell expand table
1124 ! 1999-03-01 Iredell f90 module
1125 ! 2001-02-26 Iredell ice phase
1126 !
1127 ! Usage: call gtdpi
1128 !
1129 ! Subprograms called:
1130 ! (ftdpixg) inlinable function to compute dewpoint temperature over ice
1131 !
1132 ! Attributes:
1133 ! Language: Fortran 90.
1134 !
1135 !$$$
1136  implicit none
1137  integer jx
1138  real(krealfp) xmin,xmax,xinc,t,x,pv
1139 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1140  xmin=0.1
1141  xmax=1000.1
1142  xinc=(xmax-xmin)/(nxtdpi-1)
1143  c1xtdpi=1.-xmin/xinc
1144  c2xtdpi=1./xinc
1145  t=197.0
1146  do jx=1,nxtdpi
1147  x=xmin+(jx-1)*xinc
1148  pv=x
1149  t=ftdpixg(t,pv)
1150  tbtdpi(jx)=t
1151  enddo
1152 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1153  end subroutine
1154 !-------------------------------------------------------------------------------
1155  elemental function ftdpi(pv)
1156 !$$$ Subprogram Documentation Block
1157 !
1158 ! Subprogram: ftdpi Compute dewpoint temperature over ice
1159 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1160 !
1161 ! Abstract: Compute dewpoint temperature from vapor pressure.
1162 ! A linear interpolation is done between values in a lookup table
1163 ! computed in gtdpi. See documentation for ftdpixg for details.
1164 ! Input values outside table range are reset to table extrema.
1165 ! The interpolation accuracy is better than 0.0005 Kelvin
1166 ! for dewpoint temperatures greater than 250 Kelvin,
1167 ! but decreases to 0.02 Kelvin for a dewpoint around 230 Kelvin.
1168 ! On the Cray, ftdpi is about 75 times faster than exact calculation.
1169 ! This function should be expanded inline in the calling routine.
1170 !
1171 ! Program History Log:
1172 ! 91-05-07 Iredell made into inlinable function
1173 ! 94-12-30 Iredell expand table
1174 ! 1999-03-01 Iredell f90 module
1175 ! 2001-02-26 Iredell ice phase
1176 !
1177 ! Usage: tdpi=ftdpi(pv)
1178 !
1179 ! Input argument list:
1180 ! pv Real(krealfp) vapor pressure in Pascals
1181 !
1182 ! Output argument list:
1183 ! ftdpi Real(krealfp) dewpoint temperature in Kelvin
1184 !
1185 ! Attributes:
1186 ! Language: Fortran 90.
1187 !
1188 !$$$
1189  implicit none
1190  real(krealfp) ftdpi
1191  real(krealfp),intent(in):: pv
1192  integer jx
1193  real(krealfp) xj
1194 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1195  xj=min(max(c1xtdpi+c2xtdpi*pv,1._krealfp),real(nxtdpi,krealfp))
1196  jx=min(xj,nxtdpi-1._krealfp)
1197  ftdpi=tbtdpi(jx)+(xj-jx)*(tbtdpi(jx+1)-tbtdpi(jx))
1198 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1199  end function
1200 !-------------------------------------------------------------------------------
1201  elemental function ftdpiq(pv)
1202 !$$$ Subprogram Documentation Block
1203 !
1204 ! Subprogram: ftdpiq Compute dewpoint temperature over ice
1205 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1206 !
1207 ! Abstract: Compute dewpoint temperature from vapor pressure.
1208 ! A quadratic interpolation is done between values in a lookup table
1209 ! computed in gtdpi. see documentation for ftdpixg for details.
1210 ! Input values outside table range are reset to table extrema.
1211 ! the interpolation accuracy is better than 0.00001 Kelvin
1212 ! for dewpoint temperatures greater than 250 Kelvin,
1213 ! but decreases to 0.002 Kelvin for a dewpoint around 230 Kelvin.
1214 ! On the Cray, ftdpiq is about 60 times faster than exact calculation.
1215 ! This function should be expanded inline in the calling routine.
1216 !
1217 ! Program History Log:
1218 ! 91-05-07 Iredell made into inlinable function
1219 ! 94-12-30 Iredell quadratic interpolation
1220 ! 1999-03-01 Iredell f90 module
1221 ! 2001-02-26 Iredell ice phase
1222 !
1223 ! Usage: tdpi=ftdpiq(pv)
1224 !
1225 ! Input argument list:
1226 ! pv Real(krealfp) vapor pressure in Pascals
1227 !
1228 ! Output argument list:
1229 ! ftdpiq Real(krealfp) dewpoint temperature in Kelvin
1230 !
1231 ! Attributes:
1232 ! Language: Fortran 90.
1233 !
1234 !$$$
1235  implicit none
1236  real(krealfp) ftdpiq
1237  real(krealfp),intent(in):: pv
1238  integer jx
1239  real(krealfp) xj,dxj,fj1,fj2,fj3
1240 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1241  xj=min(max(c1xtdpi+c2xtdpi*pv,1._krealfp),real(nxtdpi,krealfp))
1242  jx=min(max(nint(xj),2),nxtdpi-1)
1243  dxj=xj-jx
1244  fj1=tbtdpi(jx-1)
1245  fj2=tbtdpi(jx)
1246  fj3=tbtdpi(jx+1)
1247  ftdpiq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
1248 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1249  end function
1250 !-------------------------------------------------------------------------------
1251  elemental function ftdpix(pv)
1252 !$$$ Subprogram Documentation Block
1253 !
1254 ! Subprogram: ftdpix Compute dewpoint temperature over ice
1255 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1256 !
1257 ! Abstract: exactly compute dewpoint temperature from vapor pressure.
1258 ! An approximate dewpoint temperature for function ftdpixg
1259 ! is obtained using ftdpi so gtdpi must be already called.
1260 ! See documentation for ftdpixg for details.
1261 !
1262 ! Program History Log:
1263 ! 91-05-07 Iredell made into inlinable function
1264 ! 94-12-30 Iredell exact computation
1265 ! 1999-03-01 Iredell f90 module
1266 ! 2001-02-26 Iredell ice phase
1267 !
1268 ! Usage: tdpi=ftdpix(pv)
1269 !
1270 ! Input argument list:
1271 ! pv Real(krealfp) vapor pressure in Pascals
1272 !
1273 ! Output argument list:
1274 ! ftdpix Real(krealfp) dewpoint temperature in Kelvin
1275 !
1276 ! Subprograms called:
1277 ! (ftdpi) inlinable function to compute dewpoint temperature over ice
1278 ! (ftdpixg) inlinable function to compute dewpoint temperature over ice
1279 !
1280 ! Attributes:
1281 ! Language: Fortran 90.
1282 !
1283 !$$$
1284  implicit none
1285  real(krealfp) ftdpix
1286  real(krealfp),intent(in):: pv
1287  real(krealfp) tg
1288 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1289  tg=ftdpi(pv)
1290  ftdpix=ftdpixg(tg,pv)
1291 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1292  end function
1293 !-------------------------------------------------------------------------------
1294  elemental function ftdpixg(tg,pv)
1295 !$$$ Subprogram Documentation Block
1296 !
1297 ! Subprogram: ftdpixg Compute dewpoint temperature over ice
1298 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1299 !
1300 ! Abstract: Exactly compute dewpoint temperature from vapor pressure.
1301 ! A guess dewpoint temperature must be provided.
1302 ! The water model assumes a perfect gas, constant specific heats
1303 ! for gas and ice, and neglects the volume of the ice.
1304 ! The model does account for the variation of the latent heat
1305 ! of sublimation with temperature. The liquid option is not included.
1306 ! The Clausius-Clapeyron equation is integrated from the triple point
1307 ! to get the formula
1308 ! pvs=con_psat*(tr**xa)*exp(xb*(1.-tr))
1309 ! where tr is ttp/t and other values are physical constants.
1310 ! The formula is inverted by iterating Newtonian approximations
1311 ! for each pvs until t is found to within 1.e-6 Kelvin.
1312 ! This function can be expanded inline in the calling routine.
1313 !
1314 ! Program History Log:
1315 ! 91-05-07 Iredell made into inlinable function
1316 ! 94-12-30 Iredell exact computation
1317 ! 1999-03-01 Iredell f90 module
1318 ! 2001-02-26 Iredell ice phase
1319 !
1320 ! Usage: tdpi=ftdpixg(tg,pv)
1321 !
1322 ! Input argument list:
1323 ! tg Real(krealfp) guess dewpoint temperature in Kelvin
1324 ! pv Real(krealfp) vapor pressure in Pascals
1325 !
1326 ! Output argument list:
1327 ! ftdpixg Real(krealfp) dewpoint temperature in Kelvin
1328 !
1329 ! Attributes:
1330 ! Language: Fortran 90.
1331 !
1332 !$$$
1333  implicit none
1334  real(krealfp) ftdpixg
1335  real(krealfp),intent(in):: tg,pv
1336  real(krealfp),parameter:: terrm=1.e-6
1337  real(krealfp),parameter:: dldt=con_cvap-con_csol
1338  real(krealfp),parameter:: heat=con_hvap+con_hfus
1339  real(krealfp),parameter:: xpona=-dldt/con_rv
1340  real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
1341  real(krealfp) t,tr,pvt,el,dpvt,terr
1342  integer i
1343 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1344  t=tg
1345  do i=1,100
1346  tr=con_ttp/t
1347  pvt=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
1348  el=heat+dldt*(t-con_ttp)
1349  dpvt=el*pvt/(con_rv*t**2)
1350  terr=(pvt-pv)/dpvt
1351  t=t-terr
1352  if(abs(terr).le.terrm) exit
1353  enddo
1354  ftdpixg=t
1355 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1356  end function
1357 !-------------------------------------------------------------------------------
1358  subroutine gtdp
1359 !$$$ Subprogram Documentation Block
1360 !
1361 ! Subprogram: gtdp Compute dewpoint temperature table
1362 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1363 !
1364 ! Abstract: Compute dewpoint temperature table as a function of
1365 ! vapor pressure for inlinable function ftdp.
1366 ! Exact dewpoint temperatures are calculated in subprogram ftdpxg.
1367 ! The current implementation computes a table with a length
1368 ! of 5001 for vapor pressures ranging from 0.5 to 1000.5 Pascals
1369 ! giving a dewpoint temperature range of 208 to 319 Kelvin.
1370 !
1371 ! Program History Log:
1372 ! 91-05-07 Iredell
1373 ! 94-12-30 Iredell expand table
1374 ! 1999-03-01 Iredell f90 module
1375 ! 2001-02-26 Iredell ice phase
1376 !
1377 ! Usage: call gtdp
1378 !
1379 ! Subprograms called:
1380 ! (ftdpxg) inlinable function to compute dewpoint temperature
1381 !
1382 ! Attributes:
1383 ! Language: Fortran 90.
1384 !
1385 !$$$
1386  implicit none
1387  integer jx
1388  real(krealfp) xmin,xmax,xinc,t,x,pv
1389 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1390  xmin=0.5
1391  xmax=10000.5
1392  xinc=(xmax-xmin)/(nxtdp-1)
1393  c1xtdp=1.-xmin/xinc
1394  c2xtdp=1./xinc
1395  t=208.0
1396  do jx=1,nxtdp
1397  x=xmin+(jx-1)*xinc
1398  pv=x
1399  t=ftdpxg(t,pv)
1400  tbtdp(jx)=t
1401  enddo
1402 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1403  end subroutine
1404 !-------------------------------------------------------------------------------
1405  elemental function ftdp(pv)
1406 !$$$ Subprogram Documentation Block
1407 !
1408 ! Subprogram: ftdp Compute dewpoint temperature
1409 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1410 !
1411 ! Abstract: Compute dewpoint temperature from vapor pressure.
1412 ! A linear interpolation is done between values in a lookup table
1413 ! computed in gtdp. See documentation for ftdpxg for details.
1414 ! Input values outside table range are reset to table extrema.
1415 ! The interpolation accuracy is better than 0.0005 Kelvin
1416 ! for dewpoint temperatures greater than 250 Kelvin,
1417 ! but decreases to 0.02 Kelvin for a dewpoint around 230 Kelvin.
1418 ! On the Cray, ftdp is about 75 times faster than exact calculation.
1419 ! This function should be expanded inline in the calling routine.
1420 !
1421 ! Program History Log:
1422 ! 91-05-07 Iredell made into inlinable function
1423 ! 94-12-30 Iredell expand table
1424 ! 1999-03-01 Iredell f90 module
1425 ! 2001-02-26 Iredell ice phase
1426 !
1427 ! Usage: tdp=ftdp(pv)
1428 !
1429 ! Input argument list:
1430 ! pv Real(krealfp) vapor pressure in Pascals
1431 !
1432 ! Output argument list:
1433 ! ftdp Real(krealfp) dewpoint temperature in Kelvin
1434 !
1435 ! Attributes:
1436 ! Language: Fortran 90.
1437 !
1438 !$$$
1439  implicit none
1440  real(krealfp) ftdp
1441  real(krealfp),intent(in):: pv
1442  integer jx
1443  real(krealfp) xj
1444 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1445  xj=min(max(c1xtdp+c2xtdp*pv,1._krealfp),real(nxtdp,krealfp))
1446  jx=min(xj,nxtdp-1._krealfp)
1447  ftdp=tbtdp(jx)+(xj-jx)*(tbtdp(jx+1)-tbtdp(jx))
1448 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1449  end function
1450 !-------------------------------------------------------------------------------
1451  elemental function ftdpq(pv)
1452 !$$$ Subprogram Documentation Block
1453 !
1454 ! Subprogram: ftdpq Compute dewpoint temperature
1455 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1456 !
1457 ! Abstract: Compute dewpoint temperature from vapor pressure.
1458 ! A quadratic interpolation is done between values in a lookup table
1459 ! computed in gtdp. see documentation for ftdpxg for details.
1460 ! Input values outside table range are reset to table extrema.
1461 ! the interpolation accuracy is better than 0.00001 Kelvin
1462 ! for dewpoint temperatures greater than 250 Kelvin,
1463 ! but decreases to 0.002 Kelvin for a dewpoint around 230 Kelvin.
1464 ! On the Cray, ftdpq is about 60 times faster than exact calculation.
1465 ! This function should be expanded inline in the calling routine.
1466 !
1467 ! Program History Log:
1468 ! 91-05-07 Iredell made into inlinable function
1469 ! 94-12-30 Iredell quadratic interpolation
1470 ! 1999-03-01 Iredell f90 module
1471 ! 2001-02-26 Iredell ice phase
1472 !
1473 ! Usage: tdp=ftdpq(pv)
1474 !
1475 ! Input argument list:
1476 ! pv Real(krealfp) vapor pressure in Pascals
1477 !
1478 ! Output argument list:
1479 ! ftdpq Real(krealfp) dewpoint temperature in Kelvin
1480 !
1481 ! Attributes:
1482 ! Language: Fortran 90.
1483 !
1484 !$$$
1485  implicit none
1486  real(krealfp) ftdpq
1487  real(krealfp),intent(in):: pv
1488  integer jx
1489  real(krealfp) xj,dxj,fj1,fj2,fj3
1490 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1491  xj=min(max(c1xtdp+c2xtdp*pv,1._krealfp),real(nxtdp,krealfp))
1492  jx=min(max(nint(xj),2),nxtdp-1)
1493  dxj=xj-jx
1494  fj1=tbtdp(jx-1)
1495  fj2=tbtdp(jx)
1496  fj3=tbtdp(jx+1)
1497  ftdpq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
1498 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1499  end function
1500 !-------------------------------------------------------------------------------
1501  elemental function ftdpx(pv)
1502 !$$$ Subprogram Documentation Block
1503 !
1504 ! Subprogram: ftdpx Compute dewpoint temperature
1505 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1506 !
1507 ! Abstract: exactly compute dewpoint temperature from vapor pressure.
1508 ! An approximate dewpoint temperature for function ftdpxg
1509 ! is obtained using ftdp so gtdp must be already called.
1510 ! See documentation for ftdpxg for details.
1511 !
1512 ! Program History Log:
1513 ! 91-05-07 Iredell made into inlinable function
1514 ! 94-12-30 Iredell exact computation
1515 ! 1999-03-01 Iredell f90 module
1516 ! 2001-02-26 Iredell ice phase
1517 !
1518 ! Usage: tdp=ftdpx(pv)
1519 !
1520 ! Input argument list:
1521 ! pv Real(krealfp) vapor pressure in Pascals
1522 !
1523 ! Output argument list:
1524 ! ftdpx Real(krealfp) dewpoint temperature in Kelvin
1525 !
1526 ! Subprograms called:
1527 ! (ftdp) inlinable function to compute dewpoint temperature
1528 ! (ftdpxg) inlinable function to compute dewpoint temperature
1529 !
1530 ! Attributes:
1531 ! Language: Fortran 90.
1532 !
1533 !$$$
1534  implicit none
1535  real(krealfp) ftdpx
1536  real(krealfp),intent(in):: pv
1537  real(krealfp) tg
1538 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1539  tg=ftdp(pv)
1540  ftdpx=ftdpxg(tg,pv)
1541 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1542  end function
1543 !-------------------------------------------------------------------------------
1544  elemental function ftdpxg(tg,pv)
1545 !$$$ Subprogram Documentation Block
1546 !
1547 ! Subprogram: ftdpxg Compute dewpoint temperature
1548 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1549 !
1550 ! Abstract: Exactly compute dewpoint temperature from vapor pressure.
1551 ! A guess dewpoint temperature must be provided.
1552 ! The saturation vapor pressure over either liquid and ice is computed
1553 ! over liquid for temperatures above the triple point,
1554 ! over ice for temperatures 20 degress below the triple point,
1555 ! and a linear combination of the two for temperatures in between.
1556 ! The water model assumes a perfect gas, constant specific heats
1557 ! for gas, liquid and ice, and neglects the volume of the condensate.
1558 ! The model does account for the variation of the latent heat
1559 ! of condensation and sublimation with temperature.
1560 ! The Clausius-Clapeyron equation is integrated from the triple point
1561 ! to get the formula
1562 ! pvsl=con_psat*(tr**xa)*exp(xb*(1.-tr))
1563 ! where tr is ttp/t and other values are physical constants.
1564 ! The reference for this decision is Emanuel(1994), pages 116-117.
1565 ! The formula is inverted by iterating Newtonian approximations
1566 ! for each pvs until t is found to within 1.e-6 Kelvin.
1567 ! This function can be expanded inline in the calling routine.
1568 !
1569 ! Program History Log:
1570 ! 91-05-07 Iredell made into inlinable function
1571 ! 94-12-30 Iredell exact computation
1572 ! 1999-03-01 Iredell f90 module
1573 ! 2001-02-26 Iredell ice phase
1574 !
1575 ! Usage: tdp=ftdpxg(tg,pv)
1576 !
1577 ! Input argument list:
1578 ! tg Real(krealfp) guess dewpoint temperature in Kelvin
1579 ! pv Real(krealfp) vapor pressure in Pascals
1580 !
1581 ! Output argument list:
1582 ! ftdpxg Real(krealfp) dewpoint temperature in Kelvin
1583 !
1584 ! Attributes:
1585 ! Language: Fortran 90.
1586 !
1587 !$$$
1588  implicit none
1589  real(krealfp) ftdpxg
1590  real(krealfp),intent(in):: tg,pv
1591  real(krealfp),parameter:: terrm=1.e-6
1592  real(krealfp),parameter:: tliq=con_ttp
1593  real(krealfp),parameter:: tice=con_ttp-20.0
1594  real(krealfp),parameter:: dldtl=con_cvap-con_cliq
1595  real(krealfp),parameter:: heatl=con_hvap
1596  real(krealfp),parameter:: xponal=-dldtl/con_rv
1597  real(krealfp),parameter:: xponbl=-dldtl/con_rv+heatl/(con_rv*con_ttp)
1598  real(krealfp),parameter:: dldti=con_cvap-con_csol
1599  real(krealfp),parameter:: heati=con_hvap+con_hfus
1600  real(krealfp),parameter:: xponai=-dldti/con_rv
1601  real(krealfp),parameter:: xponbi=-dldti/con_rv+heati/(con_rv*con_ttp)
1602  real(krealfp) t,tr,w,pvtl,pvti,pvt,ell,eli,el,dpvt,terr
1603  integer i
1604 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1605  t=tg
1606  do i=1,100
1607  tr=con_ttp/t
1608  if(t.ge.tliq) then
1609  pvt=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
1610  el=heatl+dldtl*(t-con_ttp)
1611  dpvt=el*pvt/(con_rv*t**2)
1612  elseif(t.lt.tice) then
1613  pvt=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
1614  el=heati+dldti*(t-con_ttp)
1615  dpvt=el*pvt/(con_rv*t**2)
1616  else
1617  w=(t-tice)/(tliq-tice)
1618  pvtl=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
1619  pvti=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
1620  pvt=w*pvtl+(1.-w)*pvti
1621  ell=heatl+dldtl*(t-con_ttp)
1622  eli=heati+dldti*(t-con_ttp)
1623  dpvt=(w*ell*pvtl+(1.-w)*eli*pvti)/(con_rv*t**2)
1624  endif
1625  terr=(pvt-pv)/dpvt
1626  t=t-terr
1627  if(abs(terr).le.terrm) exit
1628  enddo
1629  ftdpxg=t
1630 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1631  end function
1632 !-------------------------------------------------------------------------------
1633  subroutine gthe
1634 !$$$ Subprogram Documentation Block
1635 !
1636 ! Subprogram: gthe Compute equivalent potential temperature table
1637 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1638 !
1639 ! Abstract: Compute equivalent potential temperature table
1640 ! as a function of LCL temperature and pressure over 1e5 Pa
1641 ! to the kappa power for function fthe.
1642 ! Equivalent potential temperatures are calculated in subprogram fthex
1643 ! the current implementation computes a table with a first dimension
1644 ! of 241 for temperatures ranging from 183.16 to 303.16 Kelvin
1645 ! and a second dimension of 151 for pressure over 1e5 Pa
1646 ! to the kappa power ranging from 0.04**rocp to 1.10**rocp.
1647 !
1648 ! Program History Log:
1649 ! 91-05-07 Iredell
1650 ! 94-12-30 Iredell expand table
1651 ! 1999-03-01 Iredell f90 module
1652 !
1653 ! Usage: call gthe
1654 !
1655 ! Subprograms called:
1656 ! (fthex) inlinable function to compute equiv. pot. temperature
1657 !
1658 ! Attributes:
1659 ! Language: Fortran 90.
1660 !
1661 !$$$
1662  implicit none
1663  integer jx,jy
1664  real(krealfp) xmin,xmax,ymin,ymax,xinc,yinc,x,y,pk,t
1665 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1666  xmin=con_ttp-90._krealfp
1667  xmax=con_ttp+30._krealfp
1668  ymin=0.04_krealfp**con_rocp
1669  ymax=1.10_krealfp**con_rocp
1670  xinc=(xmax-xmin)/(nxthe-1)
1671  c1xthe=1.-xmin/xinc
1672  c2xthe=1./xinc
1673  yinc=(ymax-ymin)/(nythe-1)
1674  c1ythe=1.-ymin/yinc
1675  c2ythe=1./yinc
1676  do jy=1,nythe
1677  y=ymin+(jy-1)*yinc
1678  pk=y
1679  do jx=1,nxthe
1680  x=xmin+(jx-1)*xinc
1681  t=x
1682  tbthe(jx,jy)=fthex(t,pk)
1683  enddo
1684  enddo
1685 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1686  end subroutine
1687 !-------------------------------------------------------------------------------
1688  elemental function fthe(t,pk)
1689 !$$$ Subprogram Documentation Block
1690 !
1691 ! Subprogram: fthe Compute equivalent potential temperature
1692 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1693 !
1694 ! Abstract: Compute equivalent potential temperature at the LCL
1695 ! from temperature and pressure over 1e5 Pa to the kappa power.
1696 ! A bilinear interpolation is done between values in a lookup table
1697 ! computed in gthe. see documentation for fthex for details.
1698 ! Input values outside table range are reset to table extrema,
1699 ! except zero is returned for too cold or high LCLs.
1700 ! The interpolation accuracy is better than 0.01 Kelvin.
1701 ! On the Cray, fthe is almost 6 times faster than exact calculation.
1702 ! This function should be expanded inline in the calling routine.
1703 !
1704 ! Program History Log:
1705 ! 91-05-07 Iredell made into inlinable function
1706 ! 94-12-30 Iredell expand table
1707 ! 1999-03-01 Iredell f90 module
1708 !
1709 ! Usage: the=fthe(t,pk)
1710 !
1711 ! Input argument list:
1712 ! t Real(krealfp) LCL temperature in Kelvin
1713 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
1714 !
1715 ! Output argument list:
1716 ! fthe Real(krealfp) equivalent potential temperature in Kelvin
1717 !
1718 ! Attributes:
1719 ! Language: Fortran 90.
1720 !
1721 !$$$
1722  implicit none
1723  real(krealfp) fthe
1724  real(krealfp),intent(in):: t,pk
1725  integer jx,jy
1726  real(krealfp) xj,yj,ftx1,ftx2
1727 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1728  xj=min(c1xthe+c2xthe*t,real(nxthe,krealfp))
1729  yj=min(c1ythe+c2ythe*pk,real(nythe,krealfp))
1730  if(xj.ge.1..and.yj.ge.1.) then
1731  jx=min(xj,nxthe-1._krealfp)
1732  jy=min(yj,nythe-1._krealfp)
1733  ftx1=tbthe(jx,jy)+(xj-jx)*(tbthe(jx+1,jy)-tbthe(jx,jy))
1734  ftx2=tbthe(jx,jy+1)+(xj-jx)*(tbthe(jx+1,jy+1)-tbthe(jx,jy+1))
1735  fthe=ftx1+(yj-jy)*(ftx2-ftx1)
1736  else
1737  fthe=0.
1738  endif
1739 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1740  end function
1741 !-------------------------------------------------------------------------------
1742  elemental function ftheq(t,pk)
1743 !$$$ Subprogram Documentation Block
1744 !
1745 ! Subprogram: ftheq Compute equivalent potential temperature
1746 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1747 !
1748 ! Abstract: Compute equivalent potential temperature at the LCL
1749 ! from temperature and pressure over 1e5 Pa to the kappa power.
1750 ! A biquadratic interpolation is done between values in a lookup table
1751 ! computed in gthe. see documentation for fthex for details.
1752 ! Input values outside table range are reset to table extrema,
1753 ! except zero is returned for too cold or high LCLs.
1754 ! The interpolation accuracy is better than 0.0002 Kelvin.
1755 ! On the Cray, ftheq is almost 3 times faster than exact calculation.
1756 ! This function should be expanded inline in the calling routine.
1757 !
1758 ! Program History Log:
1759 ! 91-05-07 Iredell made into inlinable function
1760 ! 94-12-30 Iredell quadratic interpolation
1761 ! 1999-03-01 Iredell f90 module
1762 !
1763 ! Usage: the=ftheq(t,pk)
1764 !
1765 ! Input argument list:
1766 ! t Real(krealfp) LCL temperature in Kelvin
1767 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
1768 !
1769 ! Output argument list:
1770 ! ftheq Real(krealfp) equivalent potential temperature in Kelvin
1771 !
1772 ! Attributes:
1773 ! Language: Fortran 90.
1774 !
1775 !$$$
1776  implicit none
1777  real(krealfp) ftheq
1778  real(krealfp),intent(in):: t,pk
1779  integer jx,jy
1780  real(krealfp) xj,yj,dxj,dyj
1781  real(krealfp) ft11,ft12,ft13,ft21,ft22,ft23,ft31,ft32,ft33
1782  real(krealfp) ftx1,ftx2,ftx3
1783 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1784  xj=min(c1xthe+c2xthe*t,real(nxthe,krealfp))
1785  yj=min(c1ythe+c2ythe*pk,real(nythe,krealfp))
1786  if(xj.ge.1..and.yj.ge.1.) then
1787  jx=min(max(nint(xj),2),nxthe-1)
1788  jy=min(max(nint(yj),2),nythe-1)
1789  dxj=xj-jx
1790  dyj=yj-jy
1791  ft11=tbthe(jx-1,jy-1)
1792  ft12=tbthe(jx-1,jy)
1793  ft13=tbthe(jx-1,jy+1)
1794  ft21=tbthe(jx,jy-1)
1795  ft22=tbthe(jx,jy)
1796  ft23=tbthe(jx,jy+1)
1797  ft31=tbthe(jx+1,jy-1)
1798  ft32=tbthe(jx+1,jy)
1799  ft33=tbthe(jx+1,jy+1)
1800  ftx1=(((ft31+ft11)/2-ft21)*dxj+(ft31-ft11)/2)*dxj+ft21
1801  ftx2=(((ft32+ft12)/2-ft22)*dxj+(ft32-ft12)/2)*dxj+ft22
1802  ftx3=(((ft33+ft13)/2-ft23)*dxj+(ft33-ft13)/2)*dxj+ft23
1803  ftheq=(((ftx3+ftx1)/2-ftx2)*dyj+(ftx3-ftx1)/2)*dyj+ftx2
1804  else
1805  ftheq=0.
1806  endif
1807 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1808  end function
1809 !-------------------------------------------------------------------------------
1810 ! elemental function fthex(t,pk)
1811  function fthex(t,pk)
1812 !$$$ Subprogram Documentation Block
1813 !
1814 ! Subprogram: fthex Compute equivalent potential temperature
1815 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1816 !
1817 ! Abstract: Exactly compute equivalent potential temperature at the LCL
1818 ! from temperature and pressure over 1e5 Pa to the kappa power.
1819 ! Equivalent potential temperature is constant for a saturated parcel
1820 ! rising adiabatically up a moist adiabat when the heat and mass
1821 ! of the condensed water are neglected. Ice is also neglected.
1822 ! The formula for equivalent potential temperature (Holton) is
1823 ! the=t*(pd**(-rocp))*exp(el*eps*pv/(cp*t*pd))
1824 ! where t is the temperature, pv is the saturated vapor pressure,
1825 ! pd is the dry pressure p-pv, el is the temperature dependent
1826 ! latent heat of condensation hvap+dldt*(t-ttp), and other values
1827 ! are physical constants defined in parameter statements in the code.
1828 ! Zero is returned if the input values make saturation impossible.
1829 ! This function should be expanded inline in the calling routine.
1830 !
1831 ! Program History Log:
1832 ! 91-05-07 Iredell made into inlinable function
1833 ! 94-12-30 Iredell exact computation
1834 ! 1999-03-01 Iredell f90 module
1835 !
1836 ! Usage: the=fthex(t,pk)
1837 !
1838 ! Input argument list:
1839 ! t Real(krealfp) LCL temperature in Kelvin
1840 ! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
1841 !
1842 ! Output argument list:
1843 ! fthex Real(krealfp) equivalent potential temperature in Kelvin
1844 !
1845 ! Attributes:
1846 ! Language: Fortran 90.
1847 !
1848 !$$$
1849  implicit none
1850  real(krealfp) fthex
1851  real(krealfp),intent(in):: t,pk
1852  real(krealfp) p,tr,pv,pd,el,expo,expmax
1853 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1854  p=pk**con_cpor
1855  tr=con_ttp/t
1856  pv=psatb*(tr**con_xpona)*exp(con_xponb*(1.-tr))
1857  pd=p-pv
1858  if(pd.gt.pv) then
1859  el=con_hvap+con_dldt*(t-con_ttp)
1860  expo=el*con_eps*pv/(con_cp*t*pd)
1861  fthex=t*pd**(-con_rocp)*exp(expo)
1862  else
1863  fthex=0.
1864  endif
1865 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1866  end function
1867 !-------------------------------------------------------------------------------
1868  subroutine gtma
1869 !$$$ Subprogram Documentation Block
1870 !
1871 ! Subprogram: gtma Compute moist adiabat tables
1872 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1873 !
1874 ! Abstract: Compute temperature and specific humidity tables
1875 ! as a function of equivalent potential temperature and
1876 ! pressure over 1e5 Pa to the kappa power for subprogram stma.
1877 ! Exact parcel temperatures are calculated in subprogram stmaxg.
1878 ! The current implementation computes a table with a first dimension
1879 ! of 151 for equivalent potential temperatures ranging from 200 to 500
1880 ! Kelvin and a second dimension of 121 for pressure over 1e5 Pa
1881 ! to the kappa power ranging from 0.01**rocp to 1.10**rocp.
1882 !
1883 ! Program History Log:
1884 ! 91-05-07 Iredell
1885 ! 94-12-30 Iredell expand table
1886 ! 1999-03-01 Iredell f90 module
1887 !
1888 ! Usage: call gtma
1889 !
1890 ! Subprograms called:
1891 ! (stmaxg) inlinable subprogram to compute parcel temperature
1892 !
1893 ! Attributes:
1894 ! Language: Fortran 90.
1895 !
1896 !$$$
1897  implicit none
1898  integer jx,jy
1899  real(krealfp) xmin,xmax,ymin,ymax,xinc,yinc,x,y,pk,the,t,q,tg
1900 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1901  xmin=200._krealfp
1902  xmax=500._krealfp
1903  ymin=0.01_krealfp**con_rocp
1904  ymax=1.10_krealfp**con_rocp
1905  xinc=(xmax-xmin)/(nxma-1)
1906  c1xma=1.-xmin/xinc
1907  c2xma=1./xinc
1908  yinc=(ymax-ymin)/(nyma-1)
1909  c1yma=1.-ymin/yinc
1910  c2yma=1./yinc
1911  do jy=1,nyma
1912  y=ymin+(jy-1)*yinc
1913  pk=y
1914  tg=xmin*y
1915  do jx=1,nxma
1916  x=xmin+(jx-1)*xinc
1917  the=x
1918  call stmaxg(tg,the,pk,t,q)
1919  tbtma(jx,jy)=t
1920  tbqma(jx,jy)=q
1921  tg=t
1922  enddo
1923  enddo
1924 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1925  end subroutine
1926 !-------------------------------------------------------------------------------
1927  elemental subroutine stma(the,pk,tma,qma)
1928 !$$$ Subprogram Documentation Block
1929 !
1930 ! Subprogram: stma Compute moist adiabat temperature
1931 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1932 !
1933 ! Abstract: Compute temperature and specific humidity of a parcel
1934 ! lifted up a moist adiabat from equivalent potential temperature
1935 ! at the LCL and pressure over 1e5 Pa to the kappa power.
1936 ! Bilinear interpolations are done between values in a lookup table
1937 ! computed in gtma. See documentation for stmaxg for details.
1938 ! Input values outside table range are reset to table extrema.
1939 ! The interpolation accuracy is better than 0.01 Kelvin
1940 ! and 5.e-6 kg/kg for temperature and humidity, respectively.
1941 ! On the Cray, stma is about 35 times faster than exact calculation.
1942 ! This subprogram should be expanded inline in the calling routine.
1943 !
1944 ! Program History Log:
1945 ! 91-05-07 Iredell made into inlinable function
1946 ! 94-12-30 Iredell expand table
1947 ! 1999-03-01 Iredell f90 module
1948 !
1949 ! Usage: call stma(the,pk,tma,qma)
1950 !
1951 ! Input argument list:
1952 ! the Real(krealfp) equivalent potential temperature in Kelvin
1953 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
1954 !
1955 ! Output argument list:
1956 ! tma Real(krealfp) parcel temperature in Kelvin
1957 ! qma Real(krealfp) parcel specific humidity in kg/kg
1958 !
1959 ! Attributes:
1960 ! Language: Fortran 90.
1961 !
1962 !$$$
1963  implicit none
1964  real(krealfp),intent(in):: the,pk
1965  real(krealfp),intent(out):: tma,qma
1966  integer jx,jy
1967  real(krealfp) xj,yj,ftx1,ftx2,qx1,qx2
1968 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1969  xj=min(max(c1xma+c2xma*the,1._krealfp),real(nxma,krealfp))
1970  yj=min(max(c1yma+c2yma*pk,1._krealfp),real(nyma,krealfp))
1971  jx=min(xj,nxma-1._krealfp)
1972  jy=min(yj,nyma-1._krealfp)
1973  ftx1=tbtma(jx,jy)+(xj-jx)*(tbtma(jx+1,jy)-tbtma(jx,jy))
1974  ftx2=tbtma(jx,jy+1)+(xj-jx)*(tbtma(jx+1,jy+1)-tbtma(jx,jy+1))
1975  tma=ftx1+(yj-jy)*(ftx2-ftx1)
1976  qx1=tbqma(jx,jy)+(xj-jx)*(tbqma(jx+1,jy)-tbqma(jx,jy))
1977  qx2=tbqma(jx,jy+1)+(xj-jx)*(tbqma(jx+1,jy+1)-tbqma(jx,jy+1))
1978  qma=qx1+(yj-jy)*(qx2-qx1)
1979 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1980  end subroutine
1981 !-------------------------------------------------------------------------------
1982  elemental subroutine stmaq(the,pk,tma,qma)
1983 !$$$ Subprogram Documentation Block
1984 !
1985 ! Subprogram: stmaq Compute moist adiabat temperature
1986 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
1987 !
1988 ! Abstract: Compute temperature and specific humidity of a parcel
1989 ! lifted up a moist adiabat from equivalent potential temperature
1990 ! at the LCL and pressure over 1e5 Pa to the kappa power.
1991 ! Biquadratic interpolations are done between values in a lookup table
1992 ! computed in gtma. See documentation for stmaxg for details.
1993 ! Input values outside table range are reset to table extrema.
1994 ! the interpolation accuracy is better than 0.0005 Kelvin
1995 ! and 1.e-7 kg/kg for temperature and humidity, respectively.
1996 ! On the Cray, stmaq is about 25 times faster than exact calculation.
1997 ! This subprogram should be expanded inline in the calling routine.
1998 !
1999 ! Program History Log:
2000 ! 91-05-07 Iredell made into inlinable function
2001 ! 94-12-30 Iredell quadratic interpolation
2002 ! 1999-03-01 Iredell f90 module
2003 !
2004 ! Usage: call stmaq(the,pk,tma,qma)
2005 !
2006 ! Input argument list:
2007 ! the Real(krealfp) equivalent potential temperature in Kelvin
2008 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
2009 !
2010 ! Output argument list:
2011 ! tmaq Real(krealfp) parcel temperature in Kelvin
2012 ! qma Real(krealfp) parcel specific humidity in kg/kg
2013 !
2014 ! Attributes:
2015 ! Language: Fortran 90.
2016 !
2017 !$$$
2018  implicit none
2019  real(krealfp),intent(in):: the,pk
2020  real(krealfp),intent(out):: tma,qma
2021  integer jx,jy
2022  real(krealfp) xj,yj,dxj,dyj
2023  real(krealfp) ft11,ft12,ft13,ft21,ft22,ft23,ft31,ft32,ft33
2024  real(krealfp) ftx1,ftx2,ftx3
2025  real(krealfp) q11,q12,q13,q21,q22,q23,q31,q32,q33,qx1,qx2,qx3
2026 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2027  xj=min(max(c1xma+c2xma*the,1._krealfp),real(nxma,krealfp))
2028  yj=min(max(c1yma+c2yma*pk,1._krealfp),real(nyma,krealfp))
2029  jx=min(max(nint(xj),2),nxma-1)
2030  jy=min(max(nint(yj),2),nyma-1)
2031  dxj=xj-jx
2032  dyj=yj-jy
2033  ft11=tbtma(jx-1,jy-1)
2034  ft12=tbtma(jx-1,jy)
2035  ft13=tbtma(jx-1,jy+1)
2036  ft21=tbtma(jx,jy-1)
2037  ft22=tbtma(jx,jy)
2038  ft23=tbtma(jx,jy+1)
2039  ft31=tbtma(jx+1,jy-1)
2040  ft32=tbtma(jx+1,jy)
2041  ft33=tbtma(jx+1,jy+1)
2042  ftx1=(((ft31+ft11)/2-ft21)*dxj+(ft31-ft11)/2)*dxj+ft21
2043  ftx2=(((ft32+ft12)/2-ft22)*dxj+(ft32-ft12)/2)*dxj+ft22
2044  ftx3=(((ft33+ft13)/2-ft23)*dxj+(ft33-ft13)/2)*dxj+ft23
2045  tma=(((ftx3+ftx1)/2-ftx2)*dyj+(ftx3-ftx1)/2)*dyj+ftx2
2046  q11=tbqma(jx-1,jy-1)
2047  q12=tbqma(jx-1,jy)
2048  q13=tbqma(jx-1,jy+1)
2049  q21=tbqma(jx,jy-1)
2050  q22=tbqma(jx,jy)
2051  q23=tbqma(jx,jy+1)
2052  q31=tbqma(jx+1,jy-1)
2053  q32=tbqma(jx+1,jy)
2054  q33=tbqma(jx+1,jy+1)
2055  qx1=(((q31+q11)/2-q21)*dxj+(q31-q11)/2)*dxj+q21
2056  qx2=(((q32+q12)/2-q22)*dxj+(q32-q12)/2)*dxj+q22
2057  qx3=(((q33+q13)/2-q23)*dxj+(q33-q13)/2)*dxj+q23
2058  qma=(((qx3+qx1)/2-qx2)*dyj+(qx3-qx1)/2)*dyj+qx2
2059 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2060  end subroutine
2061 !-------------------------------------------------------------------------------
2062  elemental subroutine stmax(the,pk,tma,qma)
2063 !$$$ Subprogram Documentation Block
2064 !
2065 ! Subprogram: stmax Compute moist adiabat temperature
2066 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2067 !
2068 ! Abstract: Exactly compute temperature and humidity of a parcel
2069 ! lifted up a moist adiabat from equivalent potential temperature
2070 ! at the LCL and pressure over 1e5 Pa to the kappa power.
2071 ! An approximate parcel temperature for subprogram stmaxg
2072 ! is obtained using stma so gtma must be already called.
2073 ! See documentation for stmaxg for details.
2074 !
2075 ! Program History Log:
2076 ! 91-05-07 Iredell made into inlinable function
2077 ! 94-12-30 Iredell exact computation
2078 ! 1999-03-01 Iredell f90 module
2079 !
2080 ! Usage: call stmax(the,pk,tma,qma)
2081 !
2082 ! Input argument list:
2083 ! the Real(krealfp) equivalent potential temperature in Kelvin
2084 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
2085 !
2086 ! Output argument list:
2087 ! tma Real(krealfp) parcel temperature in Kelvin
2088 ! qma Real(krealfp) parcel specific humidity in kg/kg
2089 !
2090 ! Subprograms called:
2091 ! (stma) inlinable subprogram to compute parcel temperature
2092 ! (stmaxg) inlinable subprogram to compute parcel temperature
2093 !
2094 ! Attributes:
2095 ! Language: Fortran 90.
2096 !
2097 !$$$
2098  implicit none
2099  real(krealfp),intent(in):: the,pk
2100  real(krealfp),intent(out):: tma,qma
2101  real(krealfp) tg,qg
2102 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2103  call stma(the,pk,tg,qg)
2104  call stmaxg(tg,the,pk,tma,qma)
2105 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2106  end subroutine
2107 !-------------------------------------------------------------------------------
2108  elemental subroutine stmaxg(tg,the,pk,tma,qma)
2109 !$$$ Subprogram Documentation Block
2110 !
2111 ! Subprogram: stmaxg Compute moist adiabat temperature
2112 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2113 !
2114 ! Abstract: exactly compute temperature and humidity of a parcel
2115 ! lifted up a moist adiabat from equivalent potential temperature
2116 ! at the LCL and pressure over 1e5 Pa to the kappa power.
2117 ! A guess parcel temperature must be provided.
2118 ! Equivalent potential temperature is constant for a saturated parcel
2119 ! rising adiabatically up a moist adiabat when the heat and mass
2120 ! of the condensed water are neglected. Ice is also neglected.
2121 ! The formula for equivalent potential temperature (Holton) is
2122 ! the=t*(pd**(-rocp))*exp(el*eps*pv/(cp*t*pd))
2123 ! where t is the temperature, pv is the saturated vapor pressure,
2124 ! pd is the dry pressure p-pv, el is the temperature dependent
2125 ! latent heat of condensation hvap+dldt*(t-ttp), and other values
2126 ! are physical constants defined in parameter statements in the code.
2127 ! The formula is inverted by iterating Newtonian approximations
2128 ! for each the and p until t is found to within 1.e-4 Kelvin.
2129 ! The specific humidity is then computed from pv and pd.
2130 ! This subprogram can be expanded inline in the calling routine.
2131 !
2132 ! Program History Log:
2133 ! 91-05-07 Iredell made into inlinable function
2134 ! 94-12-30 Iredell exact computation
2135 ! 1999-03-01 Iredell f90 module
2136 !
2137 ! Usage: call stmaxg(tg,the,pk,tma,qma)
2138 !
2139 ! Input argument list:
2140 ! tg Real(krealfp) guess parcel temperature in Kelvin
2141 ! the Real(krealfp) equivalent potential temperature in Kelvin
2142 ! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
2143 !
2144 ! Output argument list:
2145 ! tma Real(krealfp) parcel temperature in Kelvin
2146 ! qma Real(krealfp) parcel specific humidity in kg/kg
2147 !
2148 ! Attributes:
2149 ! Language: Fortran 90.
2150 !
2151 !$$$
2152  implicit none
2153  real(krealfp),intent(in):: tg,the,pk
2154  real(krealfp),intent(out):: tma,qma
2155  real(krealfp),parameter:: terrm=1.e-4
2156  real(krealfp) t,p,tr,pv,pd,el,expo,thet,dthet,terr
2157  integer i
2158 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2159  t=tg
2160  p=pk**con_cpor
2161  do i=1,100
2162  tr=con_ttp/t
2163  pv=psatb*(tr**con_xpona)*exp(con_xponb*(1.-tr))
2164  pd=p-pv
2165  el=con_hvap+con_dldt*(t-con_ttp)
2166  expo=el*con_eps*pv/(con_cp*t*pd)
2167  thet=t*pd**(-con_rocp)*exp(expo)
2168  dthet=thet/t*(1.+expo*(con_dldt*t/el+el*p/(con_rv*t*pd)))
2169  terr=(thet-the)/dthet
2170  t=t-terr
2171  if(abs(terr).le.terrm) exit
2172  enddo
2173  tma=t
2174  tr=con_ttp/t
2175  pv=psatb*(tr**con_xpona)*exp(con_xponb*(1.-tr))
2176  pd=p-pv
2177  qma=con_eps*pv/(pd+con_eps*pv)
2178 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2179  end subroutine
2180 !-------------------------------------------------------------------------------
2181  subroutine gpkap
2182 !$$$ Subprogram documentation block
2183 !
2184 ! Subprogram: gpkap Compute coefficients for p**kappa
2185 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2186 !
2187 ! Abstract: Computes pressure to the kappa table as a function of pressure
2188 ! for the table lookup function fpkap.
2189 ! Exact pressure to the kappa values are calculated in subprogram fpkapx.
2190 ! The current implementation computes a table with a length
2191 ! of 5501 for pressures ranging up to 110000 Pascals.
2192 !
2193 ! Program History Log:
2194 ! 94-12-30 Iredell
2195 ! 1999-03-01 Iredell f90 module
2196 ! 1999-03-24 Iredell table lookup
2197 !
2198 ! Usage: call gpkap
2199 !
2200 ! Subprograms called:
2201 ! fpkapx function to compute exact pressure to the kappa
2202 !
2203 ! Attributes:
2204 ! Language: Fortran 90.
2205 !
2206 !$$$
2207  implicit none
2208  integer jx
2209  real(krealfp) xmin,xmax,xinc,x,p
2210 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2211  xmin=0._krealfp
2212  xmax=110000._krealfp
2213  xinc=(xmax-xmin)/(nxpkap-1)
2214  c1xpkap=1.-xmin/xinc
2215  c2xpkap=1./xinc
2216  do jx=1,nxpkap
2217  x=xmin+(jx-1)*xinc
2218  p=x
2219  tbpkap(jx)=fpkapx(p)
2220  enddo
2221 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2222  end subroutine
2223 !-------------------------------------------------------------------------------
2224  elemental function fpkap(p)
2225 !$$$ Subprogram Documentation Block
2226 !
2227 ! Subprogram: fpkap raise pressure to the kappa power.
2228 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2229 !
2230 ! Abstract: Raise pressure over 1e5 Pa to the kappa power.
2231 ! A linear interpolation is done between values in a lookup table
2232 ! computed in gpkap. See documentation for fpkapx for details.
2233 ! Input values outside table range are reset to table extrema.
2234 ! The interpolation accuracy ranges from 9 decimal places
2235 ! at 100000 Pascals to 5 decimal places at 1000 Pascals.
2236 ! On the Cray, fpkap is over 5 times faster than exact calculation.
2237 ! This function should be expanded inline in the calling routine.
2238 !
2239 ! Program History Log:
2240 ! 91-05-07 Iredell made into inlinable function
2241 ! 94-12-30 Iredell standardized kappa,
2242 ! increased range and accuracy
2243 ! 1999-03-01 Iredell f90 module
2244 ! 1999-03-24 Iredell table lookup
2245 !
2246 ! Usage: pkap=fpkap(p)
2247 !
2248 ! Input argument list:
2249 ! p Real(krealfp) pressure in Pascals
2250 !
2251 ! Output argument list:
2252 ! fpkap Real(krealfp) p over 1e5 Pa to the kappa power
2253 !
2254 ! Attributes:
2255 ! Language: Fortran 90.
2256 !
2257 !$$$
2258  implicit none
2259  real(krealfp) fpkap
2260  real(krealfp),intent(in):: p
2261  integer jx
2262  real(krealfp) xj
2263 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2264  xj=min(max(c1xpkap+c2xpkap*p,1._krealfp),real(nxpkap,krealfp))
2265  jx=min(xj,nxpkap-1._krealfp)
2266  fpkap=tbpkap(jx)+(xj-jx)*(tbpkap(jx+1)-tbpkap(jx))
2267 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2268  end function
2269 !-------------------------------------------------------------------------------
2270  elemental function fpkapq(p)
2271 !$$$ Subprogram Documentation Block
2272 !
2273 ! Subprogram: fpkapq raise pressure to the kappa power.
2274 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2275 !
2276 ! Abstract: Raise pressure over 1e5 Pa to the kappa power.
2277 ! A quadratic interpolation is done between values in a lookup table
2278 ! computed in gpkap. see documentation for fpkapx for details.
2279 ! Input values outside table range are reset to table extrema.
2280 ! The interpolation accuracy ranges from 12 decimal places
2281 ! at 100000 Pascals to 7 decimal places at 1000 Pascals.
2282 ! On the Cray, fpkap is over 4 times faster than exact calculation.
2283 ! This function should be expanded inline in the calling routine.
2284 !
2285 ! Program History Log:
2286 ! 91-05-07 Iredell made into inlinable function
2287 ! 94-12-30 Iredell standardized kappa,
2288 ! increased range and accuracy
2289 ! 1999-03-01 Iredell f90 module
2290 ! 1999-03-24 Iredell table lookup
2291 !
2292 ! Usage: pkap=fpkapq(p)
2293 !
2294 ! Input argument list:
2295 ! p Real(krealfp) pressure in Pascals
2296 !
2297 ! Output argument list:
2298 ! fpkapq Real(krealfp) p over 1e5 Pa to the kappa power
2299 !
2300 ! Attributes:
2301 ! Language: Fortran 90.
2302 !
2303 !$$$
2304  implicit none
2305  real(krealfp) fpkapq
2306  real(krealfp),intent(in):: p
2307  integer jx
2308  real(krealfp) xj,dxj,fj1,fj2,fj3
2309 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2310  xj=min(max(c1xpkap+c2xpkap*p,1._krealfp),real(nxpkap,krealfp))
2311  jx=min(max(nint(xj),2),nxpkap-1)
2312  dxj=xj-jx
2313  fj1=tbpkap(jx-1)
2314  fj2=tbpkap(jx)
2315  fj3=tbpkap(jx+1)
2316  fpkapq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
2317 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2318  end function
2319 !-------------------------------------------------------------------------------
2320  function fpkapo(p)
2321 !$$$ Subprogram documentation block
2322 !
2323 ! Subprogram: fpkapo raise surface pressure to the kappa power.
2324 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2325 !
2326 ! Abstract: Raise surface pressure over 1e5 Pa to the kappa power
2327 ! using a rational weighted chebyshev approximation.
2328 ! The numerator is of order 2 and the denominator is of order 4.
2329 ! The pressure range is 40000-110000 Pa and kappa is defined in fpkapx.
2330 ! The accuracy of this approximation is almost 8 decimal places.
2331 ! On the Cray, fpkap is over 10 times faster than exact calculation.
2332 ! This function should be expanded inline in the calling routine.
2333 !
2334 ! Program History Log:
2335 ! 91-05-07 Iredell made into inlinable function
2336 ! 94-12-30 Iredell standardized kappa,
2337 ! increased range and accuracy
2338 ! 1999-03-01 Iredell f90 module
2339 !
2340 ! Usage: pkap=fpkapo(p)
2341 !
2342 ! Input argument list:
2343 ! p Real(krealfp) surface pressure in Pascals
2344 ! p should be in the range 40000 to 110000
2345 !
2346 ! Output argument list:
2347 ! fpkapo Real(krealfp) p over 1e5 Pa to the kappa power
2348 !
2349 ! Attributes:
2350 ! Language: Fortran 90.
2351 !
2352 !$$$
2353  implicit none
2354  real(krealfp) fpkapo
2355  real(krealfp),intent(in):: p
2356  integer,parameter:: nnpk=2,ndpk=4
2357  real(krealfp):: cnpk(0:nnpk)=(/3.13198449e-1,5.78544829e-2,&
2358  8.35491871e-4/)
2359  real(krealfp):: cdpk(0:ndpk)=(/1.,8.15968401e-2,5.72839518e-4,&
2360  -4.86959812e-7,5.24459889e-10/)
2361  integer n
2362  real(krealfp) pkpa,fnpk,fdpk
2363 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2364  pkpa=p*1.e-3_krealfp
2365  fnpk=cnpk(nnpk)
2366  do n=nnpk-1,0,-1
2367  fnpk=pkpa*fnpk+cnpk(n)
2368  enddo
2369  fdpk=cdpk(ndpk)
2370  do n=ndpk-1,0,-1
2371  fdpk=pkpa*fdpk+cdpk(n)
2372  enddo
2373  fpkapo=fnpk/fdpk
2374 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2375  end function
2376 !-------------------------------------------------------------------------------
2377  elemental function fpkapx(p)
2378 !$$$ Subprogram documentation block
2379 !
2380 ! Subprogram: fpkapx raise pressure to the kappa power.
2381 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2382 !
2383 ! Abstract: raise pressure over 1e5 Pa to the kappa power.
2384 ! Kappa is equal to rd/cp where rd and cp are physical constants.
2385 ! This function should be expanded inline in the calling routine.
2386 !
2387 ! Program History Log:
2388 ! 94-12-30 Iredell made into inlinable function
2389 ! 1999-03-01 Iredell f90 module
2390 !
2391 ! Usage: pkap=fpkapx(p)
2392 !
2393 ! Input argument list:
2394 ! p Real(krealfp) pressure in Pascals
2395 !
2396 ! Output argument list:
2397 ! fpkapx Real(krealfp) p over 1e5 Pa to the kappa power
2398 !
2399 ! Attributes:
2400 ! Language: Fortran 90.
2401 !
2402 !$$$
2403  implicit none
2404  real(krealfp) fpkapx
2405  real(krealfp),intent(in):: p
2406 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2407  fpkapx=(p/1.e5_krealfp)**con_rocp
2408 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2409  end function
2410 !-------------------------------------------------------------------------------
2411  subroutine grkap
2412 !$$$ Subprogram documentation block
2413 !
2414 ! Subprogram: grkap Compute coefficients for p**(1/kappa)
2415 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2416 !
2417 ! Abstract: Computes pressure to the 1/kappa table as a function of pressure
2418 ! for the table lookup function frkap.
2419 ! Exact pressure to the 1/kappa values are calculated in subprogram frkapx.
2420 ! The current implementation computes a table with a length
2421 ! of 5501 for pressures ranging up to 110000 Pascals.
2422 !
2423 ! Program History Log:
2424 ! 94-12-30 Iredell
2425 ! 1999-03-01 Iredell f90 module
2426 ! 1999-03-24 Iredell table lookup
2427 !
2428 ! Usage: call grkap
2429 !
2430 ! Subprograms called:
2431 ! frkapx function to compute exact pressure to the 1/kappa
2432 !
2433 ! Attributes:
2434 ! Language: Fortran 90.
2435 !
2436 !$$$
2437  implicit none
2438  integer jx
2439  real(krealfp) xmin,xmax,xinc,x,p
2440 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2441  xmin=0._krealfp
2442  xmax=fpkapx(110000._krealfp)
2443  xinc=(xmax-xmin)/(nxrkap-1)
2444  c1xrkap=1.-xmin/xinc
2445  c2xrkap=1./xinc
2446  do jx=1,nxrkap
2447  x=xmin+(jx-1)*xinc
2448  p=x
2449  tbrkap(jx)=frkapx(p)
2450  enddo
2451 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2452  end subroutine
2453 !-------------------------------------------------------------------------------
2454  elemental function frkap(pkap)
2455 !$$$ Subprogram Documentation Block
2456 !
2457 ! Subprogram: frkap raise pressure to the 1/kappa power.
2458 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2459 !
2460 ! Abstract: Raise pressure over 1e5 Pa to the 1/kappa power.
2461 ! A linear interpolation is done between values in a lookup table
2462 ! computed in grkap. See documentation for frkapx for details.
2463 ! Input values outside table range are reset to table extrema.
2464 ! The interpolation accuracy is better than 7 decimal places.
2465 ! On the IBM, fpkap is about 4 times faster than exact calculation.
2466 ! This function should be expanded inline in the calling routine.
2467 !
2468 ! Program History Log:
2469 ! 91-05-07 Iredell made into inlinable function
2470 ! 94-12-30 Iredell standardized kappa,
2471 ! increased range and accuracy
2472 ! 1999-03-01 Iredell f90 module
2473 ! 1999-03-24 Iredell table lookup
2474 !
2475 ! Usage: p=frkap(pkap)
2476 !
2477 ! Input argument list:
2478 ! pkap Real(krealfp) p over 1e5 Pa to the kappa power
2479 !
2480 ! Output argument list:
2481 ! frkap Real(krealfp) pressure in Pascals
2482 !
2483 ! Attributes:
2484 ! Language: Fortran 90.
2485 !
2486 !$$$
2487  implicit none
2488  real(krealfp) frkap
2489  real(krealfp),intent(in):: pkap
2490  integer jx
2491  real(krealfp) xj
2492 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2493  xj=min(max(c1xrkap+c2xrkap*pkap,1._krealfp),real(nxrkap,krealfp))
2494  jx=min(xj,nxrkap-1._krealfp)
2495  frkap=tbrkap(jx)+(xj-jx)*(tbrkap(jx+1)-tbrkap(jx))
2496 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2497  end function
2498 !-------------------------------------------------------------------------------
2499  elemental function frkapq(pkap)
2500 !$$$ Subprogram Documentation Block
2501 !
2502 ! Subprogram: frkapq raise pressure to the 1/kappa power.
2503 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2504 !
2505 ! Abstract: Raise pressure over 1e5 Pa to the 1/kappa power.
2506 ! A quadratic interpolation is done between values in a lookup table
2507 ! computed in grkap. see documentation for frkapx for details.
2508 ! Input values outside table range are reset to table extrema.
2509 ! The interpolation accuracy is better than 11 decimal places.
2510 ! On the IBM, fpkap is almost 4 times faster than exact calculation.
2511 ! This function should be expanded inline in the calling routine.
2512 !
2513 ! Program History Log:
2514 ! 91-05-07 Iredell made into inlinable function
2515 ! 94-12-30 Iredell standardized kappa,
2516 ! increased range and accuracy
2517 ! 1999-03-01 Iredell f90 module
2518 ! 1999-03-24 Iredell table lookup
2519 !
2520 ! Usage: p=frkapq(pkap)
2521 !
2522 ! Input argument list:
2523 ! pkap Real(krealfp) p over 1e5 Pa to the kappa power
2524 !
2525 ! Output argument list:
2526 ! frkapq Real(krealfp) pressure in Pascals
2527 !
2528 ! Attributes:
2529 ! Language: Fortran 90.
2530 !
2531 !$$$
2532  implicit none
2533  real(krealfp) frkapq
2534  real(krealfp),intent(in):: pkap
2535  integer jx
2536  real(krealfp) xj,dxj,fj1,fj2,fj3
2537 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2538  xj=min(max(c1xrkap+c2xrkap*pkap,1._krealfp),real(nxrkap,krealfp))
2539  jx=min(max(nint(xj),2),nxrkap-1)
2540  dxj=xj-jx
2541  fj1=tbrkap(jx-1)
2542  fj2=tbrkap(jx)
2543  fj3=tbrkap(jx+1)
2544  frkapq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
2545 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2546  end function
2547 !-------------------------------------------------------------------------------
2548  elemental function frkapx(pkap)
2549 !$$$ Subprogram documentation block
2550 !
2551 ! Subprogram: frkapx raise pressure to the 1/kappa power.
2552 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2553 !
2554 ! Abstract: raise pressure over 1e5 Pa to the 1/kappa power.
2555 ! Kappa is equal to rd/cp where rd and cp are physical constants.
2556 ! This function should be expanded inline in the calling routine.
2557 !
2558 ! Program History Log:
2559 ! 94-12-30 Iredell made into inlinable function
2560 ! 1999-03-01 Iredell f90 module
2561 !
2562 ! Usage: p=frkapx(pkap)
2563 !
2564 ! Input argument list:
2565 ! pkap Real(krealfp) p over 1e5 Pa to the kappa power
2566 !
2567 ! Output argument list:
2568 ! frkapx Real(krealfp) pressure in Pascals
2569 !
2570 ! Attributes:
2571 ! Language: Fortran 90.
2572 !
2573 !$$$
2574  implicit none
2575  real(krealfp) frkapx
2576  real(krealfp),intent(in):: pkap
2577 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2578  frkapx=pkap**(1/con_rocp)*1.e5_krealfp
2579 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2580  end function
2581 !-------------------------------------------------------------------------------
2582  subroutine gtlcl
2583 !$$$ Subprogram Documentation Block
2584 !
2585 ! Subprogram: gtlcl Compute equivalent potential temperature table
2586 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2587 !
2588 ! Abstract: Compute lifting condensation level temperature table
2589 ! as a function of temperature and dewpoint depression for function ftlcl.
2590 ! Lifting condensation level temperature is calculated in subprogram ftlclx
2591 ! The current implementation computes a table with a first dimension
2592 ! of 151 for temperatures ranging from 180.0 to 330.0 Kelvin
2593 ! and a second dimension of 61 for dewpoint depression ranging from
2594 ! 0 to 60 Kelvin.
2595 !
2596 ! Program History Log:
2597 ! 1999-03-01 Iredell f90 module
2598 !
2599 ! Usage: call gtlcl
2600 !
2601 ! Subprograms called:
2602 ! (ftlclx) inlinable function to compute LCL temperature
2603 !
2604 ! Attributes:
2605 ! Language: Fortran 90.
2606 !
2607 !$$$
2608  implicit none
2609  integer jx,jy
2610  real(krealfp) xmin,xmax,ymin,ymax,xinc,yinc,x,y,tdpd,t
2611 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2612  xmin=180._krealfp
2613  xmax=330._krealfp
2614  ymin=0._krealfp
2615  ymax=60._krealfp
2616  xinc=(xmax-xmin)/(nxtlcl-1)
2617  c1xtlcl=1.-xmin/xinc
2618  c2xtlcl=1./xinc
2619  yinc=(ymax-ymin)/(nytlcl-1)
2620  c1ytlcl=1.-ymin/yinc
2621  c2ytlcl=1./yinc
2622  do jy=1,nytlcl
2623  y=ymin+(jy-1)*yinc
2624  tdpd=y
2625  do jx=1,nxtlcl
2626  x=xmin+(jx-1)*xinc
2627  t=x
2628  tbtlcl(jx,jy)=ftlclx(t,tdpd)
2629  enddo
2630  enddo
2631 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2632  end subroutine
2633 !-------------------------------------------------------------------------------
2634  elemental function ftlcl(t,tdpd)
2635 !$$$ Subprogram Documentation Block
2636 !
2637 ! Subprogram: ftlcl Compute LCL temperature
2638 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2639 !
2640 ! Abstract: Compute temperature at the lifting condensation level
2641 ! from temperature and dewpoint depression.
2642 ! A bilinear interpolation is done between values in a lookup table
2643 ! computed in gtlcl. See documentation for ftlclx for details.
2644 ! Input values outside table range are reset to table extrema.
2645 ! The interpolation accuracy is better than 0.0005 Kelvin.
2646 ! On the Cray, ftlcl is ? times faster than exact calculation.
2647 ! This function should be expanded inline in the calling routine.
2648 !
2649 ! Program History Log:
2650 ! 1999-03-01 Iredell f90 module
2651 !
2652 ! Usage: tlcl=ftlcl(t,tdpd)
2653 !
2654 ! Input argument list:
2655 ! t Real(krealfp) LCL temperature in Kelvin
2656 ! tdpd Real(krealfp) dewpoint depression in Kelvin
2657 !
2658 ! Output argument list:
2659 ! ftlcl Real(krealfp) temperature at the LCL in Kelvin
2660 !
2661 ! Attributes:
2662 ! Language: Fortran 90.
2663 !
2664 !$$$
2665  implicit none
2666  real(krealfp) ftlcl
2667  real(krealfp),intent(in):: t,tdpd
2668  integer jx,jy
2669  real(krealfp) xj,yj,ftx1,ftx2
2670 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2671  xj=min(max(c1xtlcl+c2xtlcl*t,1._krealfp),real(nxtlcl,krealfp))
2672  yj=min(max(c1ytlcl+c2ytlcl*tdpd,1._krealfp),real(nytlcl,krealfp))
2673  jx=min(xj,nxtlcl-1._krealfp)
2674  jy=min(yj,nytlcl-1._krealfp)
2675  ftx1=tbtlcl(jx,jy)+(xj-jx)*(tbtlcl(jx+1,jy)-tbtlcl(jx,jy))
2676  ftx2=tbtlcl(jx,jy+1)+(xj-jx)*(tbtlcl(jx+1,jy+1)-tbtlcl(jx,jy+1))
2677  ftlcl=ftx1+(yj-jy)*(ftx2-ftx1)
2678 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2679  end function
2680 !-------------------------------------------------------------------------------
2681  elemental function ftlclq(t,tdpd)
2682 !$$$ Subprogram Documentation Block
2683 !
2684 ! Subprogram: ftlclq Compute LCL temperature
2685 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2686 !
2687 ! Abstract: Compute temperature at the lifting condensation level
2688 ! from temperature and dewpoint depression.
2689 ! A biquadratic interpolation is done between values in a lookup table
2690 ! computed in gtlcl. see documentation for ftlclx for details.
2691 ! Input values outside table range are reset to table extrema.
2692 ! The interpolation accuracy is better than 0.000003 Kelvin.
2693 ! On the Cray, ftlclq is ? times faster than exact calculation.
2694 ! This function should be expanded inline in the calling routine.
2695 !
2696 ! Program History Log:
2697 ! 1999-03-01 Iredell f90 module
2698 !
2699 ! Usage: tlcl=ftlclq(t,tdpd)
2700 !
2701 ! Input argument list:
2702 ! t Real(krealfp) LCL temperature in Kelvin
2703 ! tdpd Real(krealfp) dewpoint depression in Kelvin
2704 !
2705 ! Output argument list:
2706 ! ftlcl Real(krealfp) temperature at the LCL in Kelvin
2707 !
2708 ! Attributes:
2709 ! Language: Fortran 90.
2710 !
2711 !$$$
2712  implicit none
2713  real(krealfp) ftlclq
2714  real(krealfp),intent(in):: t,tdpd
2715  integer jx,jy
2716  real(krealfp) xj,yj,dxj,dyj
2717  real(krealfp) ft11,ft12,ft13,ft21,ft22,ft23,ft31,ft32,ft33
2718  real(krealfp) ftx1,ftx2,ftx3
2719 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2720  xj=min(max(c1xtlcl+c2xtlcl*t,1._krealfp),real(nxtlcl,krealfp))
2721  yj=min(max(c1ytlcl+c2ytlcl*tdpd,1._krealfp),real(nytlcl,krealfp))
2722  jx=min(max(nint(xj),2),nxtlcl-1)
2723  jy=min(max(nint(yj),2),nytlcl-1)
2724  dxj=xj-jx
2725  dyj=yj-jy
2726  ft11=tbtlcl(jx-1,jy-1)
2727  ft12=tbtlcl(jx-1,jy)
2728  ft13=tbtlcl(jx-1,jy+1)
2729  ft21=tbtlcl(jx,jy-1)
2730  ft22=tbtlcl(jx,jy)
2731  ft23=tbtlcl(jx,jy+1)
2732  ft31=tbtlcl(jx+1,jy-1)
2733  ft32=tbtlcl(jx+1,jy)
2734  ft33=tbtlcl(jx+1,jy+1)
2735  ftx1=(((ft31+ft11)/2-ft21)*dxj+(ft31-ft11)/2)*dxj+ft21
2736  ftx2=(((ft32+ft12)/2-ft22)*dxj+(ft32-ft12)/2)*dxj+ft22
2737  ftx3=(((ft33+ft13)/2-ft23)*dxj+(ft33-ft13)/2)*dxj+ft23
2738  ftlclq=(((ftx3+ftx1)/2-ftx2)*dyj+(ftx3-ftx1)/2)*dyj+ftx2
2739 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2740  end function
2741 !-------------------------------------------------------------------------------
2742  function ftlclo(t,tdpd)
2743 !$$$ Subprogram documentation block
2744 !
2745 ! Subprogram: ftlclo Compute LCL temperature.
2746 ! Author: Phillips org: w/NMC2X2 Date: 29 dec 82
2747 !
2748 ! Abstract: Compute temperature at the lifting condensation level
2749 ! from temperature and dewpoint depression. the formula used is
2750 ! a polynomial taken from Phillips mstadb routine which empirically
2751 ! approximates the original exact implicit relationship.
2752 ! (This kind of approximation is customary (inman, 1969), but
2753 ! the original source for this particular one is not yet known. -MI)
2754 ! Its accuracy is about 0.03 Kelvin for a dewpoint depression of 30.
2755 ! This function should be expanded inline in the calling routine.
2756 !
2757 ! Program History Log:
2758 ! 91-05-07 Iredell made into inlinable function
2759 ! 1999-03-01 Iredell f90 module
2760 !
2761 ! Usage: tlcl=ftlclo(t,tdpd)
2762 !
2763 ! Input argument list:
2764 ! t Real(krealfp) temperature in Kelvin
2765 ! tdpd Real(krealfp) dewpoint depression in Kelvin
2766 !
2767 ! Output argument list:
2768 ! ftlclo Real(krealfp) temperature at the LCL in Kelvin
2769 !
2770 ! Attributes:
2771 ! Language: Fortran 90.
2772 !
2773 !$$$
2774  implicit none
2775  real(krealfp) ftlclo
2776  real(krealfp),intent(in):: t,tdpd
2777  real(krealfp),parameter:: clcl1= 0.954442e+0,clcl2= 0.967772e-3,&
2778  clcl3=-0.710321e-3,clcl4=-0.270742e-5
2779 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2780  ftlclo=t-tdpd*(clcl1+clcl2*t+tdpd*(clcl3+clcl4*t))
2781 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2782  end function
2783 !-------------------------------------------------------------------------------
2784  elemental function ftlclx(t,tdpd)
2785 !$$$ Subprogram documentation block
2786 !
2787 ! Subprogram: ftlclx Compute LCL temperature.
2788 ! Author: Iredell org: w/NMC2X2 Date: 25 March 1999
2789 !
2790 ! Abstract: Compute temperature at the lifting condensation level
2791 ! from temperature and dewpoint depression. A parcel lifted
2792 ! adiabatically becomes saturated at the lifting condensation level.
2793 ! The water model assumes a perfect gas, constant specific heats
2794 ! for gas and liquid, and neglects the volume of the liquid.
2795 ! The model does account for the variation of the latent heat
2796 ! of condensation with temperature. The ice option is not included.
2797 ! The Clausius-Clapeyron equation is integrated from the triple point
2798 ! to get the formulas
2799 ! pvlcl=con_psat*(trlcl**xa)*exp(xb*(1.-trlcl))
2800 ! pvdew=con_psat*(trdew**xa)*exp(xb*(1.-trdew))
2801 ! where pvlcl is the saturated parcel vapor pressure at the LCL,
2802 ! pvdew is the unsaturated parcel vapor pressure initially,
2803 ! trlcl is ttp/tlcl and trdew is ttp/tdew. The adiabatic lifting
2804 ! of the parcel is represented by the following formula
2805 ! pvdew=pvlcl*(t/tlcl)**(1/kappa)
2806 ! This formula is inverted by iterating Newtonian approximations
2807 ! until tlcl is found to within 1.e-6 Kelvin. Note that the minimum
2808 ! returned temperature is 180 Kelvin.
2809 !
2810 ! Program History Log:
2811 ! 1999-03-25 Iredell
2812 !
2813 ! Usage: tlcl=ftlclx(t,tdpd)
2814 !
2815 ! Input argument list:
2816 ! t Real(krealfp) temperature in Kelvin
2817 ! tdpd Real(krealfp) dewpoint depression in Kelvin
2818 !
2819 ! Output argument list:
2820 ! ftlclx Real(krealfp) temperature at the LCL in Kelvin
2821 !
2822 ! Attributes:
2823 ! Language: Fortran 90.
2824 !
2825 !$$$
2826  implicit none
2827  real(krealfp) ftlclx
2828  real(krealfp),intent(in):: t,tdpd
2829  real(krealfp),parameter:: terrm=1.e-4,tlmin=180.,tlminx=tlmin-5.
2830  real(krealfp) tr,pvdew,tlcl,ta,pvlcl,el,dpvlcl,terr,terrp
2831  integer i
2832 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2833  tr=con_ttp/(t-tdpd)
2834  pvdew=con_psat*(tr**con_xpona)*exp(con_xponb*(1.-tr))
2835  tlcl=t-tdpd
2836  do i=1,100
2837  tr=con_ttp/tlcl
2838  ta=t/tlcl
2839  pvlcl=con_psat*(tr**con_xpona)*exp(con_xponb*(1.-tr))*ta**(1/con_rocp)
2840  el=con_hvap+con_dldt*(tlcl-con_ttp)
2841  dpvlcl=(el/(con_rv*t**2)+1/(con_rocp*tlcl))*pvlcl
2842  terr=(pvlcl-pvdew)/dpvlcl
2843  tlcl=tlcl-terr
2844  if(abs(terr).le.terrm.or.tlcl.lt.tlminx) exit
2845  enddo
2846  ftlclx=max(tlcl,tlmin)
2847 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2848  end function
2849 !-------------------------------------------------------------------------------
2850  subroutine gfuncphys
2851 !$$$ Subprogram Documentation Block
2852 !
2853 ! Subprogram: gfuncphys Compute all physics function tables
2854 ! Author: N Phillips w/NMC2X2 Date: 30 dec 82
2855 !
2856 ! Abstract: Compute all physics function tables. Lookup tables are
2857 ! set up for computing saturation vapor pressure, dewpoint temperature,
2858 ! equivalent potential temperature, moist adiabatic temperature and humidity,
2859 ! pressure to the kappa, and lifting condensation level temperature.
2860 !
2861 ! Program History Log:
2862 ! 1999-03-01 Iredell f90 module
2863 !
2864 ! Usage: call gfuncphys
2865 !
2866 ! Subprograms called:
2867 ! gpvsl compute saturation vapor pressure over liquid table
2868 ! gpvsi compute saturation vapor pressure over ice table
2869 ! gpvs compute saturation vapor pressure table
2870 ! gtdpl compute dewpoint temperature over liquid table
2871 ! gtdpi compute dewpoint temperature over ice table
2872 ! gtdp compute dewpoint temperature table
2873 ! gthe compute equivalent potential temperature table
2874 ! gtma compute moist adiabat tables
2875 ! gpkap compute pressure to the kappa table
2876 ! grkap compute pressure to the 1/kappa table
2877 ! gtlcl compute LCL temperature table
2878 !
2879 ! Attributes:
2880 ! Language: Fortran 90.
2881 !
2882 !$$$
2883  implicit none
2884 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2885  call gpvsl
2886  call gpvsi
2887  call gpvs
2888  call gtdpl
2889  call gtdpi
2890  call gtdp
2891  call gthe
2892  call gtma
2893  call gpkap
2894  call grkap
2895  call gtlcl
2896 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2897  end subroutine
2898 !-------------------------------------------------------------------------------
2899 end module
real(krealfp) c2xrkap
Definition: funcphys.f:270
real(krealfp) c1xtdp
Definition: funcphys.f:261
subroutine, public gpvsi
Definition: funcphys.f:472
real(krealfp) c2xtdpl
Definition: funcphys.f:257
integer, parameter nxma
Definition: funcphys.f:264
real(krealfp) c1xrkap
Definition: funcphys.f:270
elemental real(krealfp) function, public ftdpq(pv)
Definition: funcphys.f:1452
subroutine, public gpvs
Definition: funcphys.f:658
integer, parameter nxpvs
Definition: funcphys.f:254
real(krealfp) c1xtlcl
Definition: funcphys.f:272
real(krealfp) c2xthe
Definition: funcphys.f:263
elemental real(krealfp) function, public ftdpiq(pv)
Definition: funcphys.f:1202
real(krealfp) c2xtdp
Definition: funcphys.f:261
elemental real(krealfp) function, public frkapq(pkap)
Definition: funcphys.f:2500
real(krealfp) c2xpvs
Definition: funcphys.f:255
real(krealfp) function, public fpkapo(p)
Definition: funcphys.f:2321
elemental real(krealfp) function, public fpvsi(t)
Definition: funcphys.f:517
elemental real(krealfp) function, public ftlclx(t, tdpd)
Definition: funcphys.f:2785
elemental real(krealfp) function, public ftdp(pv)
Definition: funcphys.f:1406
integer, parameter nxpvsl
Definition: funcphys.f:250
real(krealfp) c1ythe
Definition: funcphys.f:263
elemental real(krealfp) function, public ftdpl(pv)
Definition: funcphys.f:910
elemental real(krealfp) function, public fpkapx(p)
Definition: funcphys.f:2378
elemental real(krealfp) function, public frkapx(pkap)
Definition: funcphys.f:2549
integer, parameter nxthe
Definition: funcphys.f:262
real(krealfp) c2xtlcl
Definition: funcphys.f:272
real(krealfp), dimension(nxtdp) tbtdp
Definition: funcphys.f:261
elemental subroutine, public stmaxg(tg, the, pk, tma, qma)
Definition: funcphys.f:2109
real(krealfp) c1xpvsl
Definition: funcphys.f:251
elemental real(krealfp) function, public fpvsq(t)
Definition: funcphys.f:747
real(krealfp) function, public fthex(t, pk)
Definition: funcphys.f:1812
real(krealfp) c1yma
Definition: funcphys.f:265
elemental real(krealfp) function, public ftdplx(pv)
Definition: funcphys.f:1004
real(krealfp) c2xpvsi
Definition: funcphys.f:253
subroutine, public gtma
Definition: funcphys.f:1869
real(krealfp) function, public ftlclo(t, tdpd)
Definition: funcphys.f:2743
subroutine, public gfuncphys
Definition: funcphys.f:2851
real(krealfp) c2xtdpi
Definition: funcphys.f:259
real(krealfp) c1xpvsi
Definition: funcphys.f:253
real(krealfp) c1xthe
Definition: funcphys.f:263
real(krealfp), dimension(nxpvs) tbpvs
Definition: funcphys.f:255
integer, parameter nytlcl
Definition: funcphys.f:271
subroutine, public gpkap
Definition: funcphys.f:2182
integer, parameter nxtdp
Definition: funcphys.f:260
real(krealfp) c2yma
Definition: funcphys.f:265
real(krealfp), dimension(nxtdpi) tbtdpi
Definition: funcphys.f:259
integer, parameter nyma
Definition: funcphys.f:264
elemental real(krealfp) function, public ftdpi(pv)
Definition: funcphys.f:1156
elemental real(krealfp) function, public fpvsix(t)
Definition: funcphys.f:609
elemental real(krealfp) function, public ftdpxg(tg, pv)
Definition: funcphys.f:1545
real(krealfp) c1ytlcl
Definition: funcphys.f:272
subroutine, public grkap
Definition: funcphys.f:2412
elemental real(krealfp) function, public fpvs(t)
Definition: funcphys.f:703
integer, parameter nxtdpi
Definition: funcphys.f:258
real(krealfp) c2xpvsl
Definition: funcphys.f:251
elemental real(krealfp) function, public fpvsiq(t)
Definition: funcphys.f:561
real(krealfp) c1xma
Definition: funcphys.f:265
real(krealfp), dimension(nxrkap) tbrkap
Definition: funcphys.f:270
real(krealfp), dimension(nxtlcl, nytlcl) tbtlcl
Definition: funcphys.f:272
elemental subroutine, public stmaq(the, pk, tma, qma)
Definition: funcphys.f:1983
integer, parameter nxpkap
Definition: funcphys.f:267
elemental real(krealfp) function, public ftdplxg(tg, pv)
Definition: funcphys.f:1046
real(krealfp), parameter psatb
Definition: funcphys.f:249
subroutine, public gthe
Definition: funcphys.f:1634
subroutine, public gtdpi
Definition: funcphys.f:1109
real(krealfp) c1xtdpl
Definition: funcphys.f:257
subroutine, public gtdpl
Definition: funcphys.f:864
elemental real(krealfp) function, public fpkap(p)
Definition: funcphys.f:2225
real(krealfp) c2xpkap
Definition: funcphys.f:268
real(krealfp) c2ythe
Definition: funcphys.f:263
real(krealfp) c1xpkap
Definition: funcphys.f:268
real(krealfp) c2ytlcl
Definition: funcphys.f:272
elemental subroutine, public stma(the, pk, tma, qma)
Definition: funcphys.f:1928
real(krealfp), dimension(nxpvsl) tbpvsl
Definition: funcphys.f:251
elemental real(krealfp) function, public frkap(pkap)
Definition: funcphys.f:2455
subroutine, public gtdp
Definition: funcphys.f:1359
integer, parameter nythe
Definition: funcphys.f:262
subroutine, public gpvsl
Definition: funcphys.f:290
elemental real(krealfp) function, public ftdpixg(tg, pv)
Definition: funcphys.f:1295
elemental real(krealfp) function, public fpvsx(t)
Definition: funcphys.f:795
integer, parameter nxtdpl
Definition: funcphys.f:256
subroutine, public gtlcl
Definition: funcphys.f:2583
real(krealfp), dimension(nxma, nyma) tbqma
Definition: funcphys.f:265
real(krealfp) c1xpvs
Definition: funcphys.f:255
elemental real(krealfp) function, public ftdplq(pv)
Definition: funcphys.f:955
real(krealfp), dimension(nxtdpl) tbtdpl
Definition: funcphys.f:257
integer, parameter nxpvsi
Definition: funcphys.f:252
real(krealfp) c2xma
Definition: funcphys.f:265
elemental real(krealfp) function, public ftlclq(t, tdpd)
Definition: funcphys.f:2682
real(krealfp), dimension(nxpkap) tbpkap
Definition: funcphys.f:268
elemental real(krealfp) function, public fpvslx(t)
Definition: funcphys.f:424
elemental real(krealfp) function, public ftdpx(pv)
Definition: funcphys.f:1502
elemental real(krealfp) function, public ftheq(t, pk)
Definition: funcphys.f:1743
elemental subroutine, public stmax(the, pk, tma, qma)
Definition: funcphys.f:2063
elemental real(krealfp) function, public ftlcl(t, tdpd)
Definition: funcphys.f:2635
real(krealfp), dimension(nxpvsi) tbpvsi
Definition: funcphys.f:253
real(krealfp) c1xtdpi
Definition: funcphys.f:259
elemental real(krealfp) function, public fthe(t, pk)
Definition: funcphys.f:1689
elemental real(krealfp) function, public fpvsl(t)
Definition: funcphys.f:334
integer, parameter nxtlcl
Definition: funcphys.f:271
real(krealfp), dimension(nxma, nyma) tbtma
Definition: funcphys.f:265
elemental real(krealfp) function, public ftdpix(pv)
Definition: funcphys.f:1252
integer, parameter nxrkap
Definition: funcphys.f:269
elemental real(krealfp) function, public fpvslq(t)
Definition: funcphys.f:377
integer, parameter, public krealfp
Definition: funcphys.f:246
elemental real(krealfp) function, public fpkapq(p)
Definition: funcphys.f:2271
real(krealfp), dimension(nxthe, nythe) tbthe
Definition: funcphys.f:263