This repository has been archived by the owner on Nov 8, 2024. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 51
/
module_gfs_funcphys.F
executable file
·2979 lines (2979 loc) · 110 KB
/
module_gfs_funcphys.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!-------------------------------------------------------------------------------
module module_gfs_funcphys
!$$$ Module Documentation Block
!
! Module: funcphys API for basic thermodynamic physics
! Author: Iredell Org: W/NX23 Date: 1999-03-01
!
! Abstract: This module provides an Application Program Interface
! for computing basic thermodynamic physics functions, in particular
! (1) saturation vapor pressure as a function of temperature,
! (2) dewpoint temperature as a function of vapor pressure,
! (3) equivalent potential temperature as a function of temperature
! and scaled pressure to the kappa power,
! (4) temperature and specific humidity along a moist adiabat
! as functions of equivalent potential temperature and
! scaled pressure to the kappa power,
! (5) scaled pressure to the kappa power as a function of pressure, and
! (6) temperature at the lifting condensation level as a function
! of temperature and dewpoint depression.
! The entry points required to set up lookup tables start with a "g".
! All the other entry points are functions starting with an "f" or
! are subroutines starting with an "s". These other functions and
! subroutines are elemental; that is, they return a scalar if they
! are passed only scalars, but they return an array if they are passed
! an array. These other functions and subroutines can be inlined, too.
!
! Program History Log:
! 1999-03-01 Mark Iredell
! 1999-10-15 Mark Iredell SI unit for pressure (Pascals)
! 2001-02-26 Mark Iredell Ice phase changes of Hong and Moorthi
!
! Public Variables:
! krealfp Integer parameter kind or length of reals (=kind_phys)
!
! Public Subprograms:
! gpvsl Compute saturation vapor pressure over liquid table
!
! fpvsl Elementally compute saturation vapor pressure over liquid
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvslq Elementally compute saturation vapor pressure over liquid
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvslx Elementally compute saturation vapor pressure over liquid
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! gpvsi Compute saturation vapor pressure over ice table
!
! fpvsi Elementally compute saturation vapor pressure over ice
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvsiq Elementally compute saturation vapor pressure over ice
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvsix Elementally compute saturation vapor pressure over ice
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! gpvs Compute saturation vapor pressure table
!
! fpvs Elementally compute saturation vapor pressure
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvsq Elementally compute saturation vapor pressure
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! fpvsx Elementally compute saturation vapor pressure
! function result Real(krealfp) saturation vapor pressure in Pascals
! t Real(krealfp) temperature in Kelvin
!
! gtdpl Compute dewpoint temperature over liquid table
!
! ftdpl Elementally compute dewpoint temperature over liquid
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdplq Elementally compute dewpoint temperature over liquid
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdplx Elementally compute dewpoint temperature over liquid
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdplxg Elementally compute dewpoint temperature over liquid
! function result Real(krealfp) dewpoint temperature in Kelvin
! t Real(krealfp) guess dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! gtdpi Compute dewpoint temperature table over ice
!
! ftdpi Elementally compute dewpoint temperature over ice
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpiq Elementally compute dewpoint temperature over ice
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpix Elementally compute dewpoint temperature over ice
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpixg Elementally compute dewpoint temperature over ice
! function result Real(krealfp) dewpoint temperature in Kelvin
! t Real(krealfp) guess dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! gtdp Compute dewpoint temperature table
!
! ftdp Elementally compute dewpoint temperature
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpq Elementally compute dewpoint temperature
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpx Elementally compute dewpoint temperature
! function result Real(krealfp) dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! ftdpxg Elementally compute dewpoint temperature
! function result Real(krealfp) dewpoint temperature in Kelvin
! t Real(krealfp) guess dewpoint temperature in Kelvin
! pv Real(krealfp) vapor pressure in Pascals
!
! gthe Compute equivalent potential temperature table
!
! fthe Elementally compute equivalent potential temperature
! function result Real(krealfp) equivalent potential temperature in Kelvin
! t Real(krealfp) LCL temperature in Kelvin
! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
!
! ftheq Elementally compute equivalent potential temperature
! function result Real(krealfp) equivalent potential temperature in Kelvin
! t Real(krealfp) LCL temperature in Kelvin
! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
!
! fthex Elementally compute equivalent potential temperature
! function result Real(krealfp) equivalent potential temperature in Kelvin
! t Real(krealfp) LCL temperature in Kelvin
! pk Real(krealfp) LCL pressure over 1e5 Pa to the kappa power
!
! gtma Compute moist adiabat tables
!
! stma Elementally compute moist adiabat temperature and moisture
! the Real(krealfp) equivalent potential temperature in Kelvin
! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
! tma Real(krealfp) parcel temperature in Kelvin
! qma Real(krealfp) parcel specific humidity in kg/kg
!
! stmaq Elementally compute moist adiabat temperature and moisture
! the Real(krealfp) equivalent potential temperature in Kelvin
! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
! tma Real(krealfp) parcel temperature in Kelvin
! qma Real(krealfp) parcel specific humidity in kg/kg
!
! stmax Elementally compute moist adiabat temperature and moisture
! the Real(krealfp) equivalent potential temperature in Kelvin
! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
! tma Real(krealfp) parcel temperature in Kelvin
! qma Real(krealfp) parcel specific humidity in kg/kg
!
! stmaxg Elementally compute moist adiabat temperature and moisture
! tg Real(krealfp) guess parcel temperature in Kelvin
! the Real(krealfp) equivalent potential temperature in Kelvin
! pk Real(krealfp) pressure over 1e5 Pa to the kappa power
! tma Real(krealfp) parcel temperature in Kelvin
! qma Real(krealfp) parcel specific humidity in kg/kg
!
! gpkap Compute pressure to the kappa table
!
! fpkap Elementally raise pressure to the kappa power.
! function result Real(krealfp) p over 1e5 Pa to the kappa power
! p Real(krealfp) pressure in Pascals
!
! fpkapq Elementally raise pressure to the kappa power.
! function result Real(krealfp) p over 1e5 Pa to the kappa power
! p Real(krealfp) pressure in Pascals
!
! fpkapo Elementally raise pressure to the kappa power.
! function result Real(krealfp) p over 1e5 Pa to the kappa power
! p Real(krealfp) surface pressure in Pascals
!
! fpkapx Elementally raise pressure to the kappa power.
! function result Real(krealfp) p over 1e5 Pa to the kappa power
! p Real(krealfp) pressure in Pascals
!
! grkap Compute pressure to the 1/kappa table
!
! frkap Elementally raise pressure to the 1/kappa power.
! function result Real(krealfp) pressure in Pascals
! pkap Real(krealfp) p over 1e5 Pa to the 1/kappa power
!
! frkapq Elementally raise pressure to the kappa power.
! function result Real(krealfp) pressure in Pascals
! pkap Real(krealfp) p over 1e5 Pa to the kappa power
!
! frkapx Elementally raise pressure to the kappa power.
! function result Real(krealfp) pressure in Pascals
! pkap Real(krealfp) p over 1e5 Pa to the kappa power
!
! gtlcl Compute LCL temperature table
!
! ftlcl Elementally compute LCL temperature.
! function result Real(krealfp) temperature at the LCL in Kelvin
! t Real(krealfp) temperature in Kelvin
! tdpd Real(krealfp) dewpoint depression in Kelvin
!
! ftlclq Elementally compute LCL temperature.
! function result Real(krealfp) temperature at the LCL in Kelvin
! t Real(krealfp) temperature in Kelvin
! tdpd Real(krealfp) dewpoint depression in Kelvin
!
! ftlclo Elementally compute LCL temperature.
! function result Real(krealfp) temperature at the LCL in Kelvin
! t Real(krealfp) temperature in Kelvin
! tdpd Real(krealfp) dewpoint depression in Kelvin
!
! ftlclx Elementally compute LCL temperature.
! function result Real(krealfp) temperature at the LCL in Kelvin
! t Real(krealfp) temperature in Kelvin
! tdpd Real(krealfp) dewpoint depression in Kelvin
!
! gfuncphys Compute all physics function tables
!
! Attributes:
! Language: Fortran 90
!
!$$$
use module_gfs_machine,only:kind_phys
use module_gfs_physcons
implicit none
private
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
! Public Variables
! integer,public,parameter:: krealfp=selected_real_kind(15,45)
integer,public,parameter:: krealfp=kind_phys
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
! Private Variables
real(krealfp),parameter:: psatb=con_psat*1.e-5
integer,parameter:: nxpvsl=7501
real(krealfp) c1xpvsl,c2xpvsl,tbpvsl(nxpvsl)
integer,parameter:: nxpvsi=7501
real(krealfp) c1xpvsi,c2xpvsi,tbpvsi(nxpvsi)
integer,parameter:: nxpvs=7501
real(krealfp) c1xpvs,c2xpvs,tbpvs(nxpvs)
integer,parameter:: nxtdpl=5001
real(krealfp) c1xtdpl,c2xtdpl,tbtdpl(nxtdpl)
integer,parameter:: nxtdpi=5001
real(krealfp) c1xtdpi,c2xtdpi,tbtdpi(nxtdpi)
integer,parameter:: nxtdp=5001
real(krealfp) c1xtdp,c2xtdp,tbtdp(nxtdp)
integer,parameter:: nxthe=241,nythe=151
real(krealfp) c1xthe,c2xthe,c1ythe,c2ythe,tbthe(nxthe,nythe)
integer,parameter:: nxma=151,nyma=121
real(krealfp) c1xma,c2xma,c1yma,c2yma,tbtma(nxma,nyma),tbqma(nxma,nyma)
! integer,parameter:: nxpkap=5501
integer,parameter:: nxpkap=11001
real(krealfp) c1xpkap,c2xpkap,tbpkap(nxpkap)
integer,parameter:: nxrkap=5501
real(krealfp) c1xrkap,c2xrkap,tbrkap(nxrkap)
integer,parameter:: nxtlcl=151,nytlcl=61
real(krealfp) c1xtlcl,c2xtlcl,c1ytlcl,c2ytlcl,tbtlcl(nxtlcl,nytlcl)
logical, private :: initialized=.false.
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
! Public Subprograms
public gpvsl,fpvsl,fpvslq,fpvslx
public gpvsi,fpvsi,fpvsiq,fpvsix
public gpvs,fpvs,fpvsq,fpvsx
public gtdpl,ftdpl,ftdplq,ftdplx,ftdplxg
public gtdpi,ftdpi,ftdpiq,ftdpix,ftdpixg
public gtdp,ftdp,ftdpq,ftdpx,ftdpxg
public gthe,fthe,ftheq,fthex
public gtma,stma,stmaq,stmax,stmaxg
public gpkap,fpkap,fpkapq,fpkapo,fpkapx
public grkap,frkap,frkapq,frkapx
public gtlcl,ftlcl,ftlclq,ftlclo,ftlclx
public gfuncphys
contains
!-------------------------------------------------------------------------------
subroutine gpvsl
!$$$ Subprogram Documentation Block
!
! Subprogram: gpvsl Compute saturation vapor pressure table over liquid
! Author: N Phillips W/NMC2X2 Date: 30 dec 82
!
! Abstract: Computes saturation vapor pressure table as a function of
! temperature for the table lookup function fpvsl.
! Exact saturation vapor pressures are calculated in subprogram fpvslx.
! The current implementation computes a table with a length
! of 7501 for temperatures ranging from 180. to 330. Kelvin.
!
! Program History Log:
! 91-05-07 Iredell
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
!
! Usage: call gpvsl
!
! Subprograms called:
! (fpvslx) inlinable function to compute saturation vapor pressure
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
integer jx
real(krealfp) xmin,xmax,xinc,x,t
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
xmin=180.0_krealfp
xmax=330.0_krealfp
xinc=(xmax-xmin)/(nxpvsl-1)
! c1xpvsl=1.-xmin/xinc
c2xpvsl=1./xinc
c1xpvsl=1.-xmin*c2xpvsl
do jx=1,nxpvsl
x=xmin+(jx-1)*xinc
t=x
tbpvsl(jx)=fpvslx(t)
enddo
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end subroutine
!-------------------------------------------------------------------------------
! elemental function fpvsl(t)
function fpvsl(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsl Compute saturation vapor pressure over liquid
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A linear interpolation is done between values in a lookup table
! computed in gpvsl. See documentation for fpvslx for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 6 decimal places.
! On the Cray, fpvsl is about 4 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
!
! Usage: pvsl=fpvsl(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsl Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsl
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvsl+c2xpvsl*t,1._krealfp),real(nxpvsl,krealfp))
jx=min(xj,nxpvsl-1._krealfp)
fpvsl=tbpvsl(jx)+(xj-jx)*(tbpvsl(jx+1)-tbpvsl(jx))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvslq(t)
function fpvslq(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvslq Compute saturation vapor pressure over liquid
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A quadratic interpolation is done between values in a lookup table
! computed in gpvsl. See documentation for fpvslx for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 9 decimal places.
! On the Cray, fpvslq is about 3 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell quadratic interpolation
! 1999-03-01 Iredell f90 module
!
! Usage: pvsl=fpvslq(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvslq Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvslq
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj,dxj,fj1,fj2,fj3
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvsl+c2xpvsl*t,1._krealfp),real(nxpvsl,krealfp))
jx=min(max(nint(xj),2),nxpvsl-1)
dxj=xj-jx
fj1=tbpvsl(jx-1)
fj2=tbpvsl(jx)
fj3=tbpvsl(jx+1)
fpvslq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvslx(t)
function fpvslx(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvslx Compute saturation vapor pressure over liquid
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Exactly compute saturation vapor pressure from temperature.
! The water model assumes a perfect gas, constant specific heats
! for gas and liquid, and neglects the volume of the liquid.
! The model does account for the variation of the latent heat
! of condensation with temperature. The ice option is not included.
! The Clausius-Clapeyron equation is integrated from the triple point
! to get the formula
! pvsl=con_psat*(tr**xa)*exp(xb*(1.-tr))
! where tr is ttp/t and other values are physical constants.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell exact computation
! 1999-03-01 Iredell f90 module
!
! Usage: pvsl=fpvslx(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvslx Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvslx
real(krealfp),intent(in):: t
real(krealfp),parameter:: dldt=con_cvap-con_cliq
real(krealfp),parameter:: heat=con_hvap
real(krealfp),parameter:: xpona=-dldt/con_rv
real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
real(krealfp) tr
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
tr=con_ttp/t
fpvslx=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
subroutine gpvsi
!$$$ Subprogram Documentation Block
!
! Subprogram: gpvsi Compute saturation vapor pressure table over ice
! Author: N Phillips W/NMC2X2 Date: 30 dec 82
!
! Abstract: Computes saturation vapor pressure table as a function of
! temperature for the table lookup function fpvsi.
! Exact saturation vapor pressures are calculated in subprogram fpvsix.
! The current implementation computes a table with a length
! of 7501 for temperatures ranging from 180. to 330. Kelvin.
!
! Program History Log:
! 91-05-07 Iredell
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: call gpvsi
!
! Subprograms called:
! (fpvsix) inlinable function to compute saturation vapor pressure
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
integer jx
real(krealfp) xmin,xmax,xinc,x,t
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xmin=180.0_krealfp
xmax=330.0_krealfp
xinc=(xmax-xmin)/(nxpvsi-1)
! c1xpvsi=1.-xmin/xinc
c2xpvsi=1./xinc
c1xpvsi=1.-xmin*c2xpvsi
do jx=1,nxpvsi
x=xmin+(jx-1)*xinc
t=x
tbpvsi(jx)=fpvsix(t)
enddo
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end subroutine
!-------------------------------------------------------------------------------
! elemental function fpvsi(t)
function fpvsi(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsi Compute saturation vapor pressure over ice
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A linear interpolation is done between values in a lookup table
! computed in gpvsi. See documentation for fpvsix for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 6 decimal places.
! On the Cray, fpvsi is about 4 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvsi=fpvsi(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsi Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsi
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvsi+c2xpvsi*t,1._krealfp),real(nxpvsi,krealfp))
jx=min(xj,nxpvsi-1._krealfp)
fpvsi=tbpvsi(jx)+(xj-jx)*(tbpvsi(jx+1)-tbpvsi(jx))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvsiq(t)
function fpvsiq(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsiq Compute saturation vapor pressure over ice
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A quadratic interpolation is done between values in a lookup table
! computed in gpvsi. See documentation for fpvsix for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 9 decimal places.
! On the Cray, fpvsiq is about 3 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell quadratic interpolation
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvsi=fpvsiq(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsiq Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsiq
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj,dxj,fj1,fj2,fj3
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvsi+c2xpvsi*t,1._krealfp),real(nxpvsi,krealfp))
jx=min(max(nint(xj),2),nxpvsi-1)
dxj=xj-jx
fj1=tbpvsi(jx-1)
fj2=tbpvsi(jx)
fj3=tbpvsi(jx+1)
fpvsiq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvsix(t)
function fpvsix(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsix Compute saturation vapor pressure over ice
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Exactly compute saturation vapor pressure from temperature.
! The water model assumes a perfect gas, constant specific heats
! for gas and ice, and neglects the volume of the ice.
! The model does account for the variation of the latent heat
! of condensation with temperature. The liquid option is not included.
! The Clausius-Clapeyron equation is integrated from the triple point
! to get the formula
! pvsi=con_psat*(tr**xa)*exp(xb*(1.-tr))
! where tr is ttp/t and other values are physical constants.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell exact computation
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvsi=fpvsix(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsix Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsix
real(krealfp),intent(in):: t
real(krealfp),parameter:: dldt=con_cvap-con_csol
real(krealfp),parameter:: heat=con_hvap+con_hfus
real(krealfp),parameter:: xpona=-dldt/con_rv
real(krealfp),parameter:: xponb=-dldt/con_rv+heat/(con_rv*con_ttp)
real(krealfp) tr
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
tr=con_ttp/t
fpvsix=con_psat*(tr**xpona)*exp(xponb*(1.-tr))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
subroutine gpvs
!$$$ Subprogram Documentation Block
!
! Subprogram: gpvs Compute saturation vapor pressure table
! Author: N Phillips W/NMC2X2 Date: 30 dec 82
!
! Abstract: Computes saturation vapor pressure table as a function of
! temperature for the table lookup function fpvs.
! Exact saturation vapor pressures are calculated in subprogram fpvsx.
! The current implementation computes a table with a length
! of 7501 for temperatures ranging from 180. to 330. Kelvin.
!
! Program History Log:
! 91-05-07 Iredell
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: call gpvs
!
! Subprograms called:
! (fpvsx) inlinable function to compute saturation vapor pressure
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
integer jx
real(krealfp) xmin,xmax,xinc,x,t
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xmin=180.0_krealfp
xmax=330.0_krealfp
xinc=(xmax-xmin)/(nxpvs-1)
! c1xpvs=1.-xmin/xinc
c2xpvs=1./xinc
c1xpvs=1.-xmin*c2xpvs
do jx=1,nxpvs
x=xmin+(jx-1)*xinc
t=x
tbpvs(jx)=fpvsx(t)
enddo
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end subroutine
!-------------------------------------------------------------------------------
! elemental function fpvs(t)
function fpvs(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvs Compute saturation vapor pressure
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A linear interpolation is done between values in a lookup table
! computed in gpvs. See documentation for fpvsx for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 6 decimal places.
! On the Cray, fpvs is about 4 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvs=fpvs(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvs Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvs
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvs+c2xpvs*t,1._krealfp),real(nxpvs,krealfp))
jx=min(xj,nxpvs-1._krealfp)
fpvs=tbpvs(jx)+(xj-jx)*(tbpvs(jx+1)-tbpvs(jx))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvsq(t)
function fpvsq(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsq Compute saturation vapor pressure
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute saturation vapor pressure from the temperature.
! A quadratic interpolation is done between values in a lookup table
! computed in gpvs. See documentation for fpvsx for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 9 decimal places.
! On the Cray, fpvsq is about 3 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell quadratic interpolation
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvs=fpvsq(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsq Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsq
real(krealfp),intent(in):: t
integer jx
real(krealfp) xj,dxj,fj1,fj2,fj3
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xpvs+c2xpvs*t,1._krealfp),real(nxpvs,krealfp))
jx=min(max(nint(xj),2),nxpvs-1)
dxj=xj-jx
fj1=tbpvs(jx-1)
fj2=tbpvs(jx)
fj3=tbpvs(jx+1)
fpvsq=(((fj3+fj1)/2-fj2)*dxj+(fj3-fj1)/2)*dxj+fj2
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function fpvsx(t)
function fpvsx(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsx Compute saturation vapor pressure
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Exactly compute saturation vapor pressure from temperature.
! The saturation vapor pressure over either liquid and ice is computed
! over liquid for temperatures above the triple point,
! over ice for temperatures 20 degress below the triple point,
! and a linear combination of the two for temperatures in between.
! The water model assumes a perfect gas, constant specific heats
! for gas, liquid and ice, and neglects the volume of the condensate.
! The model does account for the variation of the latent heat
! of condensation and sublimation with temperature.
! The Clausius-Clapeyron equation is integrated from the triple point
! to get the formula
! pvsl=con_psat*(tr**xa)*exp(xb*(1.-tr))
! where tr is ttp/t and other values are physical constants.
! The reference for this computation is Emanuel(1994), pages 116-117.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell exact computation
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvs=fpvsx(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsx Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) fpvsx
real(krealfp),intent(in):: t
real(krealfp),parameter:: tliq=con_ttp
real(krealfp),parameter:: tice=con_ttp-20.0
real(krealfp),parameter:: dldtl=con_cvap-con_cliq
real(krealfp),parameter:: heatl=con_hvap
real(krealfp),parameter:: xponal=-dldtl/con_rv
real(krealfp),parameter:: xponbl=-dldtl/con_rv+heatl/(con_rv*con_ttp)
real(krealfp),parameter:: dldti=con_cvap-con_csol
real(krealfp),parameter:: heati=con_hvap+con_hfus
real(krealfp),parameter:: xponai=-dldti/con_rv
real(krealfp),parameter:: xponbi=-dldti/con_rv+heati/(con_rv*con_ttp)
real(krealfp) tr,w,pvl,pvi
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
tr=con_ttp/t
if(t.ge.tliq) then
fpvsx=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
elseif(t.lt.tice) then
fpvsx=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
else
w=(t-tice)/(tliq-tice)
pvl=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
pvi=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
fpvsx=w*pvl+(1.-w)*pvi
endif
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
subroutine gtdpl
!$$$ Subprogram Documentation Block
!
! Subprogram: gtdpl Compute dewpoint temperature over liquid table
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute dewpoint temperature table as a function of
! vapor pressure for inlinable function ftdpl.
! Exact dewpoint temperatures are calculated in subprogram ftdplxg.
! The current implementation computes a table with a length
! of 5001 for vapor pressures ranging from 1 to 10001 Pascals
! giving a dewpoint temperature range of 208 to 319 Kelvin.
!
! Program History Log:
! 91-05-07 Iredell
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
!
! Usage: call gtdpl
!
! Subprograms called:
! (ftdplxg) inlinable function to compute dewpoint temperature over liquid
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
integer jx
real(krealfp) xmin,xmax,xinc,t,x,pv
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
xmin=1
xmax=10001
xinc=(xmax-xmin)/(nxtdpl-1)
c1xtdpl=1.-xmin/xinc
c2xtdpl=1./xinc
t=208.0
do jx=1,nxtdpl
x=xmin+(jx-1)*xinc
pv=x
t=ftdplxg(t,pv)
tbtdpl(jx)=t
enddo
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end subroutine
!-------------------------------------------------------------------------------
! elemental function ftdpl(pv)
function ftdpl(pv)
!$$$ Subprogram Documentation Block
!
! Subprogram: ftdpl Compute dewpoint temperature over liquid
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute dewpoint temperature from vapor pressure.
! A linear interpolation is done between values in a lookup table
! computed in gtdpl. See documentation for ftdplxg for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is better than 0.0005 Kelvin
! for dewpoint temperatures greater than 250 Kelvin,
! but decreases to 0.02 Kelvin for a dewpoint around 230 Kelvin.
! On the Cray, ftdpl is about 75 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
!
! Usage: tdpl=ftdpl(pv)
!
! Input argument list:
! pv Real(krealfp) vapor pressure in Pascals
!
! Output argument list:
! ftdpl Real(krealfp) dewpoint temperature in Kelvin
!
! Attributes:
! Language: Fortran 90.
!
!$$$
implicit none
real(krealfp) ftdpl
real(krealfp),intent(in):: pv
integer jx
real(krealfp) xj
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
if(.not.initialized) call gfuncphys()
xj=min(max(c1xtdpl+c2xtdpl*pv,1._krealfp),real(nxtdpl,krealfp))
jx=min(xj,nxtdpl-1._krealfp)
ftdpl=tbtdpl(jx)+(xj-jx)*(tbtdpl(jx+1)-tbtdpl(jx))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function
!-------------------------------------------------------------------------------
! elemental function ftdplq(pv)
function ftdplq(pv)
!$$$ Subprogram Documentation Block
!
! Subprogram: ftdplq Compute dewpoint temperature over liquid
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
! Abstract: Compute dewpoint temperature from vapor pressure.
! A quadratic interpolation is done between values in a lookup table
! computed in gtdpl. see documentation for ftdplxg for details.
! Input values outside table range are reset to table extrema.
! the interpolation accuracy is better than 0.00001 Kelvin
! for dewpoint temperatures greater than 250 Kelvin,
! but decreases to 0.002 Kelvin for a dewpoint around 230 Kelvin.
! On the Cray, ftdplq is about 60 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell quadratic interpolation
! 1999-03-01 Iredell f90 module
!
! Usage: tdpl=ftdplq(pv)
!