-
Notifications
You must be signed in to change notification settings - Fork 7
/
Copy pathdiagnostics.f90
1652 lines (1340 loc) · 59.4 KB
/
diagnostics.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
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 diagnostics
implicit none
save
private
public :: vint, hydro_interp, mix_ratio_to_spec_hum, temp_from_theta_p, &
surface_pres, density, sea_pres, lifted_index, mean_rh, &
unstag_hght, cape, sbcape, mlcape, lcl_lfc, storm_motion, &
storm_rel_hel, calcdbz, theta_phi, pvor
public :: RGas, Grav, Kappa, Eps
! Some meteorological constants (set to the values used by WRF)
real, parameter :: RGas = 287., Grav = 9.81, KappaR = 3.5, &
SpecHeatPresDry = KappaR * RGas, &
Kappa = RGas / SpecHeatPresDry, Gamma = .0065, &
LatHeatVap = 2.5e6, RVap = 461.6, Eps = RGas / RVap, &
EpsR = RVap / RGas, Missing = -9999., EpsRM1 = EpsR - 1
interface mix_ratio_to_spec_hum
module procedure mrtsh3d
module procedure mrtsh4d
end interface
interface temp_from_theta_p
module procedure tftp0d
module procedure tftp1d
module procedure tftp2d
module procedure tftp3d
module procedure tftp4d
end interface
interface density
module procedure dens3d
module procedure dens4d
end interface
contains ! Public Functions ==================================================
! vint takes a vertical profile (prof) and interpolates it to regular levels
! starting at first with an interval of delta using a profile of values for
! the interpolating vertical coordinate variable (vCoordVals) defined at the
! same points as prof. The interpolated output is returned as interpProf.
! interpType = 1 --> linear interpolation, interpType = 2 --> interpolate
! with respect to -ln(vertical coordinate value).
subroutine vint(prof, vCoordVals, first, delta, interpType, interpProf)
real, dimension(:), intent(in) :: prof, vCoordVals
integer, intent(in) :: first, delta, interpType
real, dimension(:), intent(out) :: interpProf
real, dimension(2) :: vCoordBounds, profBounds
real :: level
integer :: nz, np, k
nz = size(prof)
np = size(interpProf)
! Validate input
if (nz /= size(vCoordVals)) then
print *, "Error: Array extent mismatch in vint!"
stop
end if
! Initialize output
interpProf = Missing
! Interpolate
do k = 1, np
level = first + (k-1) * delta
! Skip this level if we don't have data for it
if (level < minval(vCoordVals) .or. level > maxval(vCoordVals)) cycle
call find_bounds(vCoordVals, prof, level, vCoordBounds, profBounds)
if (interpType == 2) then ! Handle logarithmic interpolation case
vCoordBounds = -log(vCoordBounds)
level = -log(level)
end if
interpProf(k) = lin_int(vCoordBounds(1), vCoordBounds(2), &
profBounds(1), profBounds(2), level)
end do
end subroutine vint
! ===========================================================================
! hydro_interp hydrostatically interpolates geopotential heights on pressure
! surfaces underground. Provided as input are the height (zs; m), pressure
! (ps; hPa), and temperature (ts; K) in the lowest WRF model layer, as well
! as the lowest (in altitude) pressure surface desired (pBot) and the
! pressure level spacing (dp). The routine takes the currently computed
! geopotential heights on pressure surfaces (z), and interpolates to those
! grid points determined to be underground.
subroutine hydro_interp(zs, ps, ts, pBot, dp, z)
real, dimension(:,:,:), intent(in) :: zs, ps, ts
integer, intent(in) :: pBot, dp
real, dimension(:,:,:,:), intent(inout) :: z
real, parameter :: LapseRate = .0065, Expo = RGas * LapseRate / Grav
integer :: nz, k, p
nz = size(z,3)
! Validate
if (any(shape(zs) /= shape(ps)) .or. any(shape(zs) /= shape(ts)) .or. &
size(z,1) /= size(zs,1) .or. size(z,2) /= size(zs,2) .or. &
size(z,4) /= size(zs,3)) stop "Array mismatch in hydro_interp!"
do k = 1, nz
if (all(z /= Missing)) exit ! Are we above the highest terrain?
p = pBot + (k-1) * dp
where (z(:,:,k,:) == Missing)
z(:,:,k,:) = zs - ts / LapseRate * ((p/ps)**Expo - 1)
end where
end do
end subroutine hydro_interp
! ===========================================================================
! Both mrtsh3d and mrtsh4d convert mixing ratio (kg/kg) to specific
! humidity (kg/kg).
function mrtsh3d(w) result (q)
real, dimension(:,:,:), intent(in) :: w
real, dimension(size(w,1),size(w,2),size(w,3)) :: q
q = w / (1 + w)
end function mrtsh3d
! ...........................................................................
function mrtsh4d(w) result (q)
real, dimension(:,:,:,:), intent(in) :: w
real, dimension(size(w,1),size(w,2),size(w,3),size(w,4)) :: q
q = w / (1 + w)
end function mrtsh4d
! ===========================================================================
! tftp0d, tftp1d, tftp2d, tftp3d, and tftp4d calculate the temperature (K)
! given the potential temperature (th; K) and pressure (p; Pa or hPa
! depending on si flag).
function tftp0d(th, p, si) result (t)
real, intent(in) :: th, p
logical, intent(in), optional :: si
real :: t
real :: p0r
! Set appropriate p0
p0r = 1e-5 ! Default is SI units (p0 = 100000)
if (present(si)) then
if (.not. si) p0r = .001 ! (p0 = 1000)
end if
! Calculate temperature (exp and log faster than ** on my machine)
t = th * exp(Kappa * log(p * p0r))
end function tftp0d
! ...........................................................................
function tftp1d(th, p, si) result (t)
real, dimension(:), intent(in) :: th, p
logical, intent(in), optional :: si
real, dimension(size(p,1)) :: t
real :: p0r
! Make sure arrays match up
if (any(shape(th) /= shape(p))) stop "Array extent mismatch in tftp2d!"
! Set appropriate p0
p0r = 1e-5 ! Default is SI units (p0 = 100000)
if (present(si)) then
if (.not. si) p0r = .001 ! (p0 = 1000)
end if
! Calculate temperature
t = th * (p * p0r) ** Kappa
end function tftp1d
! ...........................................................................
function tftp2d(th, p, si) result (t)
real, dimension(:,:), intent(in) :: th, p
logical, intent(in), optional :: si
real, dimension(size(p,1),size(p,2)) :: t
real :: p0r
! Make sure arrays match up
if (any(shape(th) /= shape(p))) stop "Array extent mismatch in tftp2d!"
! Set appropriate p0
p0r = 1e-5 ! Default is SI units (p0 = 100000)
if (present(si)) then
if (.not. si) p0r = .001 ! (p0 = 1000)
end if
! Calculate temperature
t = th * (p * p0r) ** Kappa
end function tftp2d
! ...........................................................................
function tftp3d(th, p, si) result (t)
real, dimension(:,:,:), intent(in) :: th, p
logical, intent(in), optional :: si
real, dimension(size(p,1),size(p,2),size(p,3)) :: t
real :: p0r
! Make sure arrays match up
if (any(shape(th) /= shape(p))) stop "Array extent mismatch in tftp3d!"
! Set appropriate p0
p0r = 1e-5 ! Default is SI units (p0 = 100000)
if (present(si)) then
if (.not. si) p0r = .001 ! (p0 = 1000)
end if
! Calculate temperature
t = th * (p * p0r) ** Kappa
end function tftp3d
! ...........................................................................
function tftp4d(th, p, si) result (t)
real, dimension(:,:,:,:), intent(in) :: th, p
logical, intent(in), optional :: si
real, dimension(size(p,1),size(p,2),size(p,3),size(p,4)) :: t
real :: p0r
! Make sure arrays match up
if (any(shape(th) /= shape(p))) stop "Array extent mismatch in tftp4d!"
! Set appropriate p0
p0r = 1e-5 ! Default is SI units (p0 = 100000)
if (present(si)) then
if (.not. si) p0r = .001 ! (p0 = 1000)
end if
! Calculate temperature
t = th * (p * p0r) ** Kappa
end function tftp4d
! ===========================================================================
! surface_pres calculates the surface pressure (Pa) given the following data:
! ht- Geopotential height (m) at lowest two w (full) levels (x,y,z,t)
! p1- Pressure (Pa) at lowest mass (half) level (x,y,t)
! q1- Specific humidity (kg/kg) at lowest mass level (x,y,t)
! t1- Temperature (K) at lowest mass level (x,y,t)
! sfcq- Specific humidity (kq/kq) at surface/2-m above ground (x,y,t)
! It is assumed that the array extents match appropriately.
function surface_pres(ht, p1, q1, t1, sfcq)
real, dimension(:,:,:,:), intent(in) :: ht
real, dimension(:,:,:), intent(in) :: p1, q1, t1, sfcq
real, dimension(size(p1,1),size(p1,2),size(p1,3)) :: surface_pres
real, dimension(size(t1,1),size(t1,2),size(t1,3)) :: tvirt
tvirt = (1 + .608 * .5 * (sfcq + q1)) * t1
surface_pres = p1 * exp((Grav / (RGas * tvirt)) &
* .5 * (ht(:,:,2,:) - ht(:,:,1,:)))
end function surface_pres
! ===========================================================================
! dens3d and dens4d calculate the density (kg/m**3) using the ideal gas law
! given:
! p- Pressure (Pa)
! t- Temperature (K)
! q- Water vapor mixing ratio (kg/kg)
! It is assumed that the array extents match appropriately.
function dens3d(p, t, q)
real, dimension(:,:,:), intent(in) :: p, t, q
real, dimension(size(p,1),size(p,2),size(p,3)) :: dens3d
dens3d = p / (t * (RGas + q * RVap))
dens3d = dens3d * (1 + q)
end function dens3d
! ...........................................................................
function dens4d(p, t, q)
real, dimension(:,:,:,:), intent(in) :: p, t, q
real, dimension(size(p,1),size(p,2),size(p,3),size(p,4)) :: dens4d
dens4d = p / (t * (RGas + q * RVap))
dens4d = dens4d * (1 + q)
end function dens4d
! ===========================================================================
! sea_pres is taken from the WRF NCL tools. The inputs are:
! zs- Geopotential height (m) on w (full) levels (x,y,z,t)
! theta- Potential temperature (K) on mass (half) levels (x,y,z,t)
! p- Pressure (Pa) on mass levels (x,y,z,t)
! q- Specific humidity (kg/kg) on mass levels (x,y,z,t)
! It is assumed that the array extents match appropriately.
function sea_pres(zs, theta, p, q)
! Estimate sea level pressure.
real, dimension(:,:,:,:), intent(in) :: zs, theta, p, q
real, dimension(size(p,1),size(p,2),size(p,4)) :: sea_pres
real, dimension(size(p,1),size(p,2),size(p,3),size(p,4)) :: t, z
real, dimension(size(p,1),size(p,2)) :: tSurf, tSeaLevel
integer, dimension(size(p,1),size(p,2)) :: level
! Specific constants for assumptions made in this routine:
real, parameter :: Tc = 273.16+17.5, PConst = 10000
! Local variables:
real :: plo, phi, tlo, thi, zlo, zhi, pAtPConst, tAtPConst, zAtPConst
integer :: nx, ny, nz, nt
integer :: i, j, k, klo, khi, l
logical :: l1, l2, l3, found
nx = size(p,1)
ny = size(p,2)
nz = size(p,3)
nt = size(p,4)
! Convert potential temperature to temperature
t = temp_from_theta_p(theta, p)
! Unstagger z
z(:,:,1:nz,:) = .5 * (zs(:,:,1:nz,:) + zs(:,:,2:nz+1,:))
! Grand loop over all times
do l = 1, nt
! Find least zeta level that is PConst Pa above the surface. We later
! use this level to extrapolate a surface pressure and temperature,
! which is supposed to reduce the effect of the diurnal heating cycle in
! the pressure field.
do j = 1 , ny
do i = 1 , nx
level(i,j) = -1
k = 1
found = .false.
do while( (.not. found) .and. (k <= nz))
if ( p(i,j,k,l) < p(i,j,1,l)-PConst ) then
level(i,j) = k
found = .true.
end if
k = k+1
end do
if ( level(i,j) == -1 ) then
print '(A,I4,A)','Troubles finding level ', &
nint(PConst)/100,' above ground.'
print '(A,I4,A,I4,A)', &
'Problems first occur at (',i,',',j,')'
print '(A,F6.1,A)', &
'Surface pressure = ',p(i,j,1,l)/100,' hPa.'
stop 'Error_in_finding_100_hPa_up'
end if
end do
end do
! Get temperature PConst Pa above surface. Use this to extrapolate
! the temperature at the surface and down to sea level.
do j = 1 , ny
do i = 1 , nx
klo = max ( level(i,j) - 1 , 1 )
khi = min ( klo + 1 , nz - 1 )
if ( klo == khi ) then
print '(A)','Trapping levels are weird.'
print '(A,I3,A,I3,A)','klo = ',klo,', khi = ',khi, &
': and they should not be equal.'
stop 'Error_trapping_levels'
end if
plo = p(i,j,klo,l)
phi = p(i,j,khi,l)
tlo = t(i,j,klo,l) * (1. + 0.608 * q(i,j,klo,l) )
thi = t(i,j,khi,l) * (1. + 0.608 * q(i,j,khi,l) )
zlo = z(i,j,klo,l)
zhi = z(i,j,khi,l)
pAtPConst = p(i,j,1,l) - PConst
TAtPConst = thi - (thi-tlo) * log(pAtPConst/phi) * log(plo/phi)
zAtPConst = zhi - (zhi-zlo) * log(pAtPConst/phi) * log(plo/phi)
tSurf(i,j) = tAtPConst * (p(i,j,1,l)/pAtPConst)**(Gamma*RGas/Grav)
tSeaLevel(i,j) = tAtPConst + Gamma * zAtPConst
end do
end do
! If we follow a traditional computation, there is a correction to the
! sea level temperature if both the surface and sea level temnperatures
! are *too* hot.
do j = 1 , ny
do i = 1 , nx
l1 = tSeaLevel(i,j) < Tc
l2 = tSurf (i,j) <= Tc
l3 = .not. l1
if ( l2 .and. l3 ) then
tSeaLevel(i,j) = Tc
else
tSeaLevel(i,j) = Tc - 0.005 * (tSurf(i,j) - Tc)**2
end if
end do
end do
! The grand finale: ta da!
sea_pres(:,:,l) = p(:,:,1,l) * &
exp((2.*Grav*z(:,:,1,l)) / (RGas*(tSeaLevel + tSurf)))
end do
end function sea_pres
! ===========================================================================
! lifted_index calculates the lifted index (K) based on the supplied profiles
! of potential temperature (theta; K), pressure (p; hPa), and mixing ratio
! (w; kg/kg). Array indices are ordered from the ground up. The lowest 4
! model levels are mixed to produce the hypothetical parcel, which is then
! lifted adiabatically to saturation, and pseudoadiabatically thereafter.
function lifted_index(theta, p, w) result (li)
real, dimension(:), intent(in) :: theta, p, w
real :: li
real, dimension(2) :: pBounds, thBounds
real :: dens, thMix, wMix, sm, tk, pLCL, thes, tParcel, tEnv
integer :: nz, k
nz = size(theta)
! Check that array extents match
if (nz /= size(p) .or. nz /= size(w)) &
stop "Array mismatch in lifted_index."
! Set lifted index to 0 if we don't have enough data
if (nz < 4 .or. p(1) < 500.) then
li = 0
else
! Determine temperature of pseudoadiabatically lifted parcel at 500 mb
! Step 1: Mix the lowest 4 model levels
thMix = 0; wMix = 0; sm = 0
do k = 1, 4
dens = density_thp(theta(k), p(k)) ! Effect of water vapor ignored
thMix = thMix + dens * theta(k)
wMix = wMix + dens * w(k)
sm = sm + dens
end do
thMix = thMix / sm
wMix = wMix / sm
! Step 2: Lift parcel to LCL
tk = temp_from_theta_p(thMix, p(1), si=.false.)
! This equation is derived from Curry & Webster (1999, p. 174)
pLCL = 1000 * (temp_lcl(tk, rh(wMix, p(1), tk)) / thMix)**KappaR
! Step 3: Lift parcel to 500 mb
thes = theta_es(thMix, pLCL)
tParcel = temp_from_theta_p(theta_es_to_theta_fast(thes, 500.), 500., &
si=.false.)
! Get environmental temperature at 500 mb
call find_bounds(p, theta, 500., pBounds, thBounds)
tEnv = temp_from_theta_p(lin_int(pBounds(1), pBounds(2), thBounds(1), &
thBounds(2), 500.), 500., si=.false.)
li = tEnv - tParcel
end if
end function lifted_index
! ===========================================================================
! mean_rh computes the mean relative humidity (%) between pBot hPa and pTop
! hPa based on the supplied profiles of potential temperature (theta; K),
! pressure (p; hPa), and mixing ratio (w; kg/kg). Also provided is the
! surface pressure (pSfc; hPa). Array indices are ordered from the ground
! up. The mean is weighted by the thickness of each layer.
function mean_rh(theta, p, w, pSfc, pBot, pTop) result (mrh)
real, dimension(:), intent(in) :: theta, p, w
real, intent(in) :: pSfc, pBot, pTop
real :: mrh
real, dimension(0:size(theta)-1) :: pMid
real :: pThick, relHum, totPresThick, totRH
integer :: nz, k
nz = size(theta)
! Check that array extents match
if (nz /= size(p) .or. nz /= size(w)) stop "Array mismatch in mean_rh."
! Calculate pressure at midpoints
pMid(0) = pSfc
pMid(1:nz-1) = .5 * (p(1:nz-1) + p(2:nz))
! Calculate mean RH
totPresThick = 0; totRH = 0
do k = 1, nz-1
if (pMid(k) > pBot) cycle
if (pMid(k-1) < pTop) exit
! Calculate layer thickness
if (pMid(k-1) > pBot) then
pThick = pBot - pMid(k)
else if (pMid(k) < pTop) then
pThick = pMid(k-1) - pTop
else
pThick = pMid(k-1) - pMid(k)
end if
relHum = rh(w(k), p(k), temp_from_theta_p(theta(k), p(k), si=.false.))
totRH = totRH + pThick * relHum
totPresThick = totPresThick + pThick
end do
if (totPresThick > .01) then
mrh = 100 * totRH / totPresThick
else
mrh = 0
end if
end function mean_rh
! ===========================================================================
! unstag_hght unstaggers the geopotential height from w levels (zs) given
! by znw to u levels (z) given by znu. mu and pt provide additional
! information needed to do an interpolation based on ln p.
function unstag_hght(zs, znw, znu, mu, pt) result (z)
real, dimension(:), intent(in) :: zs, znw, znu
real, intent(in) :: mu, pt
real, dimension(size(znu)) :: z
real, dimension(size(znw)) :: pw
real, dimension(size(znu)) :: pu
integer :: k
! Validate input
if (size(zs) /= size(znw) .or. size(znw) /= size(znu) + 1) &
print *, "Warning: Array extent mismatch in unstag_hght!"
! Initialize pressures
do k = 1, size(znw)
if (k < size(znw)) pu(k) = pd(znu(k))
pw(k) = pd(znw(k))
end do
! Interpolate to half levels
do k = 1, size(znu)
!orig z(k) = lin_int(-log(pw(k)), -log(pw(k+1)), zs(k), zs(k+1), -log(pu(k)))
z(k) = lin_int(pw(k), pw(k+1), zs(k), zs(k+1), pu(k))
end do
contains
real function pd(eta)
real, intent(in) :: eta
pd = pt + mu * eta
end function pd
end function unstag_hght
! ===========================================================================
! cape calculates the convective available potential energy (J/kg) given
! profiles of potential temperature (theta; K), pressure (p; hPa),
! mixing ratio (w; kg/kg), and height (z; m). By default, it returns the
! most unstable CAPE, but if lplOutI is set to .true., the height above
! ground level of the parcel whose CAPE is most unstable is returned instead.
! In that case, the terrain height at the profile location (ter) is also
! required.
real function cape(theta, p, w, z, ter, lplOutI)
real, dimension(:), intent(in) :: theta, p, w, z
real, optional, intent(in) :: ter
logical, optional, intent(in) :: lplOutI
real, dimension(size(z)) :: tProf
real, dimension(2) :: pBounds, zBounds
real :: lpl, tk, pLCL, thes, zLCL, parcelCape
integer :: nz, i
logical :: lplOut
! Initialize
cape = 0; lpl = 0
nz = size(z)
lplOut = .false.
if (present(lplOutI)) then
lplOut = lplOutI
if (.not. present(ter)) stop "For LPL, need terrain height!"
end if
! Validate
if (nz /= size(theta) .or. nz /= size(p) .or. nz /= size(w)) &
stop "Error: Array extent mismatch in cape!"
! Compute environmental profile with virtual temperature correction.
! Note that applying the correction here doesn't affect LCL values below.
tProf = temp_from_theta_p(theta, p, si=.false.) * (1 + EpsRM1 * w)
! We loop over all parcels originating in the lowest 300 hPa
do i = 1, nz-1
if (p(i) < p(1) - 300) exit
! Find LCL for this parcel
tk = temp_from_theta_p(theta(i), p(i), si=.false.)
pLCL = 1000 * (temp_lcl(tk, rh(w(i), p(i), tk)) / theta(i))**KappaR
if (pLCL > p(i)) pLCL = p(i)
thes = theta_es(theta(i), pLCL)
call find_bounds(p, z, pLCL, pBounds, zBounds)
zLCL = lin_int(pBounds(1), pBounds(2), zBounds(1), zBounds(2), pLCL)
parcelCape = calc_cape(tProf(i:nz), p(i:nz), z(i:nz), theta(i), w(i), &
thes, pLCL, zLCL)
! Check if largest CAPE so far
if (parcelCape > cape) then
cape = parcelCape
lpl = z(i) - ter
end if
end do
if (lplOut) cape = lpl
end function cape
! ===========================================================================
! sbcape calculates the convective available potential energy (J/kg) of the
! parcel lifted from the lowest model layer given profiles of potential
! temperature (theta; K), pressure (p; hPa), mixing ratio (w; kg/kg), and
! height (z; m). By default, it returns the CAPE, but if cinOutI is set to
! .true., the convective inhibition (J/kg) is returned instead.
real function sbcape(theta, p, w, z, cinOutI)
real, dimension(:), intent(in) :: theta, p, w, z
logical, optional, intent(in) :: cinOutI
real, dimension(size(z)) :: tProf
real, dimension(2) :: pBounds, zBounds
real :: tk, pLCL, zLCL, thes
integer :: nz
logical :: cinOut
! Initialize
sbcape = 0
nz = size(z)
cinOut = .false.
if (present(cinOutI)) cinOut = cinOutI
! Validate
if (nz /= size(theta) .or. nz /= size(p) .or. nz /= size(w)) &
stop "Error: Array extent mismatch in sbcape!"
! Compute environmental profile with virtual temperature correction.
! Note that applying the correction here doesn't affect LCL values below.
tProf = temp_from_theta_p(theta, p, si=.false.) * (1 + EpsRM1 * w)
! Find LCL for lowest parcel
tk = temp_from_theta_p(theta(1), p(1), si=.false.)
pLCL = 1000 * (temp_lcl(tk, rh(w(1), p(1), tk)) / theta(1))**KappaR
if (pLCL > p(1)) pLCL = p(1)
thes = theta_es(theta(1), pLCL)
! Calculate what we want
if (cinOut) then
sbcape = calc_cin(tProf, p, z, theta(1), w(1), thes, pLCL)
else
call find_bounds(p, z, pLCL, pBounds, zBounds)
zLCL = lin_int(pBounds(1), pBounds(2), zBounds(1), zBounds(2), pLCL)
sbcape = calc_cape(tProf, p, z, theta(1), w(1), thes, pLCL, zLCL)
end if
end function sbcape
! ===========================================================================
! mlcape calculates the convective available potential energy (J/kg) of the
! parcel mixed from the lowest 4 model layers given profiles of potential
! temperature (theta; K), pressure (p; hPa), mixing ratio (w; kg/kg), and
! height (z; m). By default, it returns the CAPE, but if cinOutI is set to
! .true., the convective inhibition (J/kg) is returned instead.
real function mlcape(theta, p, w, z, cinOutI)
real, dimension(:), intent(in) :: theta, p, w, z
logical, optional, intent(in) :: cinOutI
real, dimension(size(z)) :: tProf
real, dimension(2) :: pBounds, zBounds
real :: thMix, wMix, total, dens, tk, pLCL, zLCL, thes
integer :: nz, k
logical :: cinOut
! Initialize
mlcape = 0
nz = size(z)
cinOut = .false.
if (present(cinOutI)) cinOut = cinOutI
! Validate
if (nz /= size(theta) .or. nz /= size(p) .or. nz /= size(w)) &
stop "Error: Array extent mismatch in mlcape!"
! Set result to 0 if we don't have enough data
if (nz < 4) then
mlcape = 0
return
end if
! Mix the lowest 4 model levels
thMix = 0; wMix = 0; total = 0
do k = 1, 4
dens = density_thp(theta(k), p(k)) ! Effect of water vapor ignored
thMix = thMix + dens * theta(k)
wMix = wMix + dens * w(k)
total = total + dens
end do
thMix = thMix / total
wMix = wMix / total
! Compute environmental profile with virtual temperature correction.
! Note that applying the correction here doesn't affect LCL values below.
tProf = temp_from_theta_p(theta, p, si=.false.) * (1 + EpsRM1 * w)
! Find LCL for mixed parcel
tk = temp_from_theta_p(thMix, p(1), si=.false.)
pLCL = 1000 * (temp_lcl(tk, rh(wMix, p(1), tk)) / thMix)**KappaR
if (pLCL > p(1)) pLCL = p(1)
thes = theta_es(thMix, pLCL)
! Calculate what we want
if (cinOut) then
mlcape = calc_cin(tProf, p, z, thMix, wMix, thes, pLCL)
else
call find_bounds(p, z, pLCL, pBounds, zBounds)
zLCL = lin_int(pBounds(1), pBounds(2), zBounds(1), zBounds(2), pLCL)
mlcape = calc_cape(tProf, p, z, thMix, wMix, thes, pLCL, zLCL)
end if
end function mlcape
! ===========================================================================
! lcl_lfc returns the lifting condensation level or level of free convection.
! If lclFlag is true (false), the former (latter) is returned. Provided as
! input are profiles of potential temperature (theta; K), pressure (p; hPa),
! mixing ratio (w; kg/kg), and geopotential height (z; m). The terrain
! height (ter; m) is also provided as input.
real function lcl_lfc(theta, p, w, z, ter, lclFlag)
real, dimension(:), intent(in) :: theta, p, w, z
real, intent(in) :: ter
logical, intent(in) :: lclFlag
real, dimension(size(z)) :: tProf, tTraj
real, dimension(2) :: pBounds, zBounds
real :: thMix, wMix, total, dens, tk, pLCL, zLCL, thes
integer :: nz, k
! Initialize
lcl_lfc = -9999 ! The GEMPAK missing value
nz = size(z)
! Validate
if (nz /= size(theta) .or. nz /= size(p) .or. nz /= size(w)) &
stop "Error: Array extent mismatch in lcl_lfc!"
! Return immediately if we don't have enough data
if (nz < 4) return
! Mix the lowest 4 model levels
thMix = 0; wMix = 0; total = 0
do k = 1, 4
dens = density_thp(theta(k), p(k)) ! Effect of water vapor ignored
thMix = thMix + dens * theta(k)
wMix = wMix + dens * w(k)
total = total + dens
end do
thMix = thMix / total
wMix = wMix / total
! Find LCL for mixed parcel
tk = temp_from_theta_p(thMix, p(1), si=.false.)
pLCL = 1000 * (temp_lcl(tk, rh(wMix, p(1), tk)) / thMix)**KappaR
if (pLCL > p(1)) pLCL = p(1)
call find_bounds(p, z, pLCL, pBounds, zBounds)
zLCL = lin_int(pBounds(1), pBounds(2), zBounds(1), zBounds(2), pLCL)
! Calculate what we want
if (lclFlag) then
lcl_lfc = zLCL - ter
else
thes = theta_es(thMix, pLCL)
! Compute environmental profile and ascent trajectory using the
! virtual temperature correction.
tProf = temp_from_theta_p(theta, p, si=.false.) * (1 + EpsRM1 * w)
tTraj = calc_traj(p, tProf(1), theta(1), w(1), thes, pLCL)
! Find LFC height (from bottom up)
do k = 2, nz-1
if (p(k) >= pLCL .or. tTraj(k) <= tProf(k)) cycle
! We just passed through LFC, but we need to allow for case where the
! LCL is in a region of positive area already due to a superadiabatic
! environmental lapse rate.
if (tTraj(k-1) - tProf(k-1) >= 0) then
lcl_lfc = zLCL
else
lcl_lfc = lin_int(tTraj(k-1)-tProf(k-1), tTraj(k)-tProf(k), &
z(k-1), z(k), 0.)
end if
exit
end do
if (lcl_lfc /= -9999) lcl_lfc = lcl_lfc - ter
end if
end function lcl_lfc
! ===========================================================================
! storm_motion computes either the u or v component of the supercell storm
! motion vector derived by Bunkers et al. (2000). uComp = true (false)
! means the u (v) component will be output. Provided for input are the u and
! v components of the wind (m), the 10-m winds (u and v; m), and the
! terrain height (ter; m).
function storm_motion(u, v, z, u10m, v10m, ter, uComp)
real, dimension(:,:,:,:), intent(in) :: u, v, z
real, dimension(:,:,:), intent(in) :: u10m, v10m
real, dimension(:,:), intent(in) :: ter
logical, intent(in) :: uComp
real, dimension(size(u,1),size(u,2),size(u,4)) :: storm_motion
real, parameter :: D = 7.5
integer, parameter :: NumLevs = 13, Dz = 500 ! 0 to 6 km
real, dimension(NumLevs), parameter :: Heights = (/ (Dz * k, &
k = 0, NumLevs-1) /)
real, dimension(3), parameter :: HeadHeights = (/ 5500, 5750, 6000 /), &
TailHeights = (/ 0, 250, 500 /)
real, dimension(size(u,1),size(u,2),size(u,4)) :: uShear, vShear, magShear
integer :: k, nx, ny, nt
nx = size(u,1); ny = size(u,2); nt = size(u,4)
! Validate
if (any(shape(u) /= shape(z)) .or. any(shape(v) /= shape(z)) .or. &
nx /= size(u10m,1) .or. ny /= size(u10m,2) .or. &
nt /= size(u10m,3) .or. nx /= size(v10m,1) .or. &
ny /= size(v10m,2) .or. nt /= size(v10m,3) .or. &
nx /= size(ter,1) .or. ny /= size(ter,2)) &
stop "Array mismatch in storm_motion!"
! Calculate Bunkers wind shear and its magnitude
uShear = mean_wind(u, z, u10m, ter, HeadHeights) - &
mean_wind(u, z, u10m, ter, TailHeights)
vShear = mean_wind(v, z, v10m, ter, HeadHeights) - &
mean_wind(v, z, v10m, ter, TailHeights)
magShear = sqrt(uShear**2 + vShear**2)
! Calculate appropriate storm motion vector.
if (uComp) then
storm_motion = mean_wind(u, z, u10m, ter, Heights) + &
D * VShear / MagShear
else
storm_motion = mean_wind(v, z, v10m, ter, Heights) - &
D * UShear / MagShear
end if
end function storm_motion
! ===========================================================================
! storm_rel_hel computes the storm relative helicity (m2/s2) for the layer
! defined by heights (m), given as input the u and v components of the wind
! (m), the geopotential height (z; m), the u and v components of the wind at
! 10m (u10m, v10m; m), and the terrain height (ter; m).
function storm_rel_hel(u, v, z, u10m, v10m, ter, heights) result (srh)
real, dimension(:,:,:,:), intent(in) :: u, v, z
real, dimension(:,:,:), intent(in) :: u10m, v10m
real, dimension(:,:), intent(in) :: ter
real, dimension(:), intent(in) :: heights
real, dimension(size(u,1),size(u,2),size(u,4)) :: srh
real, dimension(size(z,1),size(z,2),size(z,3),size(z,4)) :: hagl
real, dimension(size(u,1),size(u,2),size(heights),size(u,4)) :: uInt, &
vInt, &
xVort, &
yVort, &
srhs
real, dimension(size(u,1),size(u,2),size(u,4)) :: uStorm, vStorm
real, dimension(2) :: haglBounds, windBounds
integer :: nx, ny, nz, nt, numLevs, i, j, k, t
nx = size(u,1); ny = size(u,2); nz = size(u,3); nt = size(u,4)
numLevs = size(heights)
! Validate
if (any(shape(u) /= shape(z)) .or. any(shape(v) /= shape(z)) .or. &
nx /= size(u10m,1) .or. ny /= size(u10m,2) .or. &
nt /= size(u10m,3) .or. nx /= size(v10m,1) .or. &
ny /= size(v10m,2) .or. nt /= size(v10m,3) .or. &
nx /= size(ter,1) .or. ny /= size(ter,2)) &
stop "Array mismatch in storm_motion!"
hagl = z - &
spread(source=spread(source=ter, dim=3, ncopies=nz), dim=4, ncopies=nt)
! Interpolate winds to constant height grid
do t = 1, nt
do j = 1, ny
do i = 1, nx
do k = 1, numLevs
if (heights(k) < hagl(i,j,1,t)) then
uInt(i,j,k,t) = u10m(i,j,t)
vInt(i,j,k,t) = v10m(i,j,t)
else
call find_bounds(hagl(i,j,:,t), u(i,j,:,t), heights(k), &
haglBounds, windBounds)
uInt(i,j,k,t) = lin_int(haglBounds(1), haglBounds(2), &
windBounds(1), windBounds(2), heights(k))
call find_bounds(hagl(i,j,:,t), v(i,j,:,t), heights(k), &
haglBounds, windBounds)
vInt(i,j,k,t) = lin_int(haglBounds(1), haglBounds(2), &
windBounds(1), windBounds(2), heights(k))
end if
end do
end do
end do
end do
! Neglect vertical motion component of horizontal vorticity, a valid
! approximation on the synoptic and mesoscales.
xVort = -calc_ddz(vInt, heights)
yVort = calc_ddz(uInt, heights)
uStorm = storm_motion(u, v, z, u10m, v10m, ter, .true.)
vStorm = storm_motion(u, v, z, u10m, v10m, ter, .false.)
srhs = (u - spread(uStorm, dim=3, ncopies=numLevs)) * xVort + &
(v - spread(vStorm, dim=3, ncopies=numLevs)) * yVort
srh = 0
do k = 1, numLevs-1
srh = srh + .5*(heights(k+1)-heights(k))*(srhs(:,:,k,:)+srhs(:,:,k+1,:))
end do
end function storm_rel_hel
! ===========================================================================
! calcdbz diagnoses reflectivity (dBZ) provided inputs of air density
! (rhoair; kg/m^3), temperature (tmk; K), rain water mixing ratio
! (qra; kg/kg), snow mixing ratio (qsn; kg/kg), and graupel mixing ratio
! (qgr; kg/kg). This function is largely taken from that in RIP4. Note that
! to be more accurate the parameters in this function should depend on the
! microphysical scheme being used.
function calcdbz(rhoair, tmk, qra, qsn, qgr) result (dbz)
real, dimension(:,:,:,:), intent(in) :: rhoair, tmk, qra, qsn, qgr
real, dimension(size(tmk,1),size(tmk,2),size(tmk,3),size(tmk,4)) :: dbz
real, parameter :: r1 = 1.e-15, ron2 = 1.e10, gon = 5.e7, ron_min = 8.e6, &
ron_qr0 = 0.00010, ron_delqr0 = 0.25*ron_qr0, &
ron_const1r = (ron2 - ron_min)*0.5, &
ron_const2r = (ron2 + ron_min)*0.5, &
gamma_seven = 720., rhowat = 1000., &
rho_r = rhowat, & ! 1000. kg m^-3
rho_s = 100., & ! kg m^-3
rho_g = 400., & ! kg m^-3
alpha = 0.224, celkel = 273.15
real, dimension(size(tmk,1),size(tmk,2),size(tmk,3),size(tmk,4)) :: ronv, &
sonv, &
gonv
real :: pi, factor_r, factor_s, factor_g
! Constants
pi = 4. * atan(1.)
factor_r = gamma_seven * 1.e18 * (1./(pi*rho_r))**1.75
factor_s = gamma_seven * 1.e18 * (1./(pi*rho_s))**1.75 &
* (rho_s/rhowat)**2 * alpha
factor_g = gamma_seven * 1.e18 * (1./(pi*rho_g))**1.75 &
* (rho_g/rhowat)**2 * alpha
! Validate
if (any(shape(tmk) /= shape(rhoair)) .or. any(shape(tmk) /= shape(qra)) &
.or. any(shape(tmk) /= shape(qsn)) &
.or. any(shape(tmk) /= shape(qgr))) stop "Array mismatch in calcdbz!"
! Calculate variable intercept parameters
sonv = min(2.0e8, 2.0e6 * exp(-0.12 * min(-0.001, tmk-celkel)))
where (qgr <= r1)
gonv = gon
elsewhere
gonv = max(1.e4, min(2.38 * (pi*rho_g / (rhoair*qgr))**0.92, gon))
end where
where (qra <= r1)
ronv = ron2
elsewhere