-
Notifications
You must be signed in to change notification settings - Fork 112
/
params.py
3383 lines (2841 loc) · 118 KB
/
params.py
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
# -*- coding: utf-8 -*-
''' Thermodynamic Parameter Routines '''
from __future__ import division
import numpy as np
import numpy.ma as ma
from sharppy.sharptab import interp, utils, thermo, winds
from sharppy.sharptab.constants import *
'''
This file contains various functions to perform the calculation of various convection indices.
Because of this, parcel lifting routines are also found in this file.
Functions denoted with a (*) in the docstring refer to functions that were added to the SHARPpy package that
were not ported from the Storm Prediction Center. They have been included as they have been used by the
community in an effort to expand SHARPpy to support the many parameters used in atmospheric science.
While the logic for these functions are based in the scientific literature, validation
of the output from these functions is occasionally difficult to perform. Although we have made an effort
to resolve code issues when they arise, values from these functions may be erronious and may require additional
inspection by the user. We appreciate any contributions by the meteorological community that can
help better validate these SHARPpy functions!
'''
__all__ = ['DefineParcel', 'Parcel', 'inferred_temp_advection']
__all__ += ['k_index', 't_totals', 'c_totals', 'v_totals', 'precip_water']
__all__ += ['temp_lvl', 'max_temp', 'mean_mixratio', 'mean_theta', 'mean_thetae', 'mean_relh']
__all__ += ['lapse_rate', 'max_lapse_rate', 'most_unstable_level', 'parcelx', 'bulk_rich']
__all__ += ['bunkers_storm_motion', 'effective_inflow_layer']
__all__ += ['convective_temp', 'esp', 'pbl_top', 'precip_eff', 'dcape', 'sig_severe']
__all__ += ['dgz', 'ship', 'stp_cin', 'stp_fixed', 'scp', 'mmp', 'wndg', 'sherb', 'tei', 'cape']
__all__ += ['mburst', 'dcp', 'ehi', 'sweat', 'hgz', 'lhp', 'integrate_parcel']
class DefineParcel(object):
'''
Create a parcel from a supplied profile object.
Parameters
----------
prof : profile object
Profile object
Optional Keywords
flag : int (default = 1)
Parcel Selection
1 - Observed Surface Parcel
2 - Forecast Surface Parcel
3 - Most Unstable Parcel
4 - Mean Mixed Layer Parcel
5 - User Defined Parcel
6 - Mean Effective Layer Parcel
Optional Keywords (Depending on Parcel Selected)
Parcel (flag) == 1: Observed Surface Parcel
None
Parcel (flag) == 2: Forecast Surface Parcel
pres : number (default = 100 hPa)
Depth over which to mix the boundary layer; only changes
temperature; does not affect moisture
Parcel (flag) == 3: Most Unstable Parcel
pres : number (default = 400 hPa)
Depth over which to look for the the most unstable parcel
starting from the surface pressure
Parcel (flag) == 4: Mixed Layer Parcel
pres : number (default = 100 hPa)
Depth over which to mix the surface parcel
Parcel (flag) == 5: User Defined Parcel
pres : number (default = SFC - 100 hPa)
Pressure of the parcel to lift
tmpc : number (default = Temperature at the provided pressure)
Temperature of the parcel to lift
dwpc : number (default = Dew Point at the provided pressure)
Dew Point of the parcel to lift
Parcel (flag) == 6: Effective Inflow Layer
ecape : number (default = 100)
The minimum amount of CAPE a parcel needs to be considered
part of the inflow layer
ecinh : number (default = -250)
The maximum amount of CINH allowed for a parcel to be
considered as part of the inflow layer
'''
def __init__(self, prof, flag, **kwargs):
self.flag = flag
if flag == 1:
self.presval = prof.pres[prof.sfc]
self.__sfc(prof)
elif flag == 2:
self.presval = kwargs.get('pres', 100)
self.__fcst(prof, **kwargs)
elif flag == 3:
self.presval = kwargs.get('pres', 300)
self.__mu(prof, **kwargs)
elif flag == 4:
self.presval = kwargs.get('pres', 100)
self.__ml(prof, **kwargs)
elif flag == 5:
self.presval = kwargs.get('pres', prof.pres[prof.sfc])
self.__user(prof, **kwargs)
elif flag == 6:
self.presval = kwargs.get('pres', 100)
self.__effective(prof, **kwargs)
else:
self.presval = kwargs.get('pres', prof.gSndg[prof.sfc])
self.__sfc(prof)
def __sfc(self, prof):
'''
Create a parcel using surface conditions
'''
self.desc = 'Surface Parcel'
self.pres = prof.pres[prof.sfc]
self.tmpc = prof.tmpc[prof.sfc]
self.dwpc = prof.dwpc[prof.sfc]
def __fcst(self, prof, **kwargs):
'''
Create a parcel using forecast conditions.
'''
self.desc = 'Forecast Surface Parcel'
self.tmpc = max_temp(prof)
self.pres = prof.pres[prof.sfc]
pbot = self.pres; ptop = self.pres - 100.
self.dwpc = thermo.temp_at_mixrat(mean_mixratio(prof, pbot, ptop, exact=True), self.pres)
def __mu(self, prof, **kwargs):
'''
Create the most unstable parcel within the lowest XXX hPa, where
XXX is supplied. Default XXX is 400 hPa.
'''
self.desc = 'Most Unstable Parcel in Lowest %.2f hPa' % self.presval
pbot = prof.pres[prof.sfc]
ptop = pbot - self.presval
self.pres = most_unstable_level(prof, pbot=pbot, ptop=ptop)
self.tmpc = interp.temp(prof, self.pres)
self.dwpc = interp.dwpt(prof, self.pres)
def __ml(self, prof, **kwargs):
'''
Create a mixed-layer parcel with mixing within the lowest XXX hPa,
where XXX is supplied. Default is 100 hPa.
If
'''
self.desc = '%.2f hPa Mixed Layer Parcel' % self.presval
pbot = kwargs.get('pbot', prof.pres[prof.sfc])
ptop = pbot - self.presval
self.pres = pbot
mtheta = mean_theta(prof, pbot, ptop, exact=True)
self.tmpc = thermo.theta(1000., mtheta, self.pres)
mmr = mean_mixratio(prof, pbot, ptop, exact=True)
self.dwpc = thermo.temp_at_mixrat(mmr, self.pres)
def __user(self, prof, **kwargs):
'''
Create a user defined parcel.
'''
self.desc = '%.2f hPa Parcel' % self.presval
self.pres = self.presval
self.tmpc = kwargs.get('tmpc', interp.temp(prof, self.pres))
self.dwpc = kwargs.get('dwpc', interp.dwpt(prof, self.pres))
def __effective(self, prof, **kwargs):
'''
Create the mean-effective layer parcel.
'''
ecape = kwargs.get('ecape', 100)
ecinh = kwargs.get('ecinh', -250)
pbot, ptop = effective_inflow_layer(prof, ecape, ecinh)
if utils.QC(pbot) and pbot > 0:
self.desc = '%.2f hPa Mean Effective Layer Centered at %.2f' % ( pbot-ptop, (pbot+ptop)/2.)
mtha = mean_theta(prof, pbot, ptop)
mmr = mean_mixratio(prof, pbot, ptop)
self.pres = (pbot+ptop)/2.
self.tmpc = thermo.theta(1000., mtha, self.pres)
self.dwpc = thermo.temp_at_mixrat(mmr, self.pres)
else:
self.desc = 'Defaulting to Surface Layer'
self.pres = prof.pres[prof.sfc]
self.tmpc = prof.tmpc[prof.sfc]
self.dwpc = prof.dwpc[prof.sfc]
if utils.QC(pbot): self.pbot = pbot
else: self.pbot = ma.masked
if utils.QC(ptop): self.ptop = ptop
else: self.pbot = ma.masked
class Parcel(object):
'''
Initialize the parcel variables
Parameters
----------
pbot : number
Lower-bound (pressure; hPa) that the parcel is lifted
ptop : number
Upper-bound (pressure; hPa) that the parcel is lifted
pres : number
Pressure of the parcel to lift (hPa)
tmpc : number
Temperature of the parcel to lift (C)
dwpc : number
Dew Point of the parcel to lift (C)
Attributes
----------
pres : number
parcel beginning pressure (mb)
tmpc : number
parcel beginning temperature (C)
dwpc : number
parcel beginning dewpoint (C)
ptrace : array
parcel trace pressure (mb)
ttrace : array
parcel trace temperature (C)
blayer : number
Pressure of the bottom of the layer the parcel is lifted (mb)
tlayer : number
Pressure of the top of the layer the parcel is lifted (mb)
entrain : number
Parcel entrainment fraction (not yet implemented)
lclpres : number
Parcel LCL (lifted condensation level) pressure (mb)
lclhght : number
Parcel LCL height (m AGL)
lfcpres : number
Parcel LFC (level of free convection) pressure (mb)
lfchght: number
Parcel LCL height (m AGL)
elpres : number
Parcel EL (equilibrium level) pressure (mb)
elhght : number
Parcel EL height (m AGL)
mplpres : number
Maximum Parcel Level (mb)
mplhght : number
Maximum Parcel Level (m AGL)
bplus : number
Parcel CAPE (J/kg)
bminus : number
Parcel CIN below 500 mb (J/kg)
bfzl : number
Parcel CAPE up to freezing level (J/kg)
b3km : number
Parcel CAPE up to 3 km (J/kg)
b6km : number
Parcel CAPE up to 6 km (J/kg)
p0c: number
Pressure value at 0 C (mb)
pm10c : number
Pressure value at -10 C (mb)
pm20c : number
Pressure value at -20 C (mb)
pm30c : number
Pressure value at -30 C (mb)
hght0c : number
Height value at 0 C (m AGL)
hghtm10c : number
Height value at -10 C (m AGL)
hghtm20c : number
Height value at -20 C (m AGL)
hghtm30c : number
Height value at -30 C (m AGL)
wm10c : number
Wetbulb at -10 C (C)
wm20c : number
Wetbulb at -20 C (C)
wm30c : number
Wetbulb at -30 C (C)
li5 : number
500-mb lifted index (C)
li3 : number
300-mb lifted index (C)
brnshear : number
Bulk Richardson Number Shear (kts)
brnu : number
U-component Bulk Richardson Number Shear (kts)
brnv : number
V-component Bulk Richardson Number Shear (kts)
brn : number
Bulk Richardson Number (unitless)
limax : number
Maximum lifted index value (C)
limaxpres : number
Pressure at Maximum lifted index (mb)
cap : number
Cap strength (C)
cappres : number
Cap strength pressure (mb)
bmin : number
Buoyancy minimum (C)
bminpres : number
Pressure at the buoyancy minimum (mb)
'''
def __init__(self, **kwargs):
self.pres = ma.masked # Parcel beginning pressure (mb)
self.tmpc = ma.masked # Parcel beginning temperature (C)
self.dwpc = ma.masked # Parcel beginning dewpoint (C)
self.ptrace = ma.masked # Parcel trace pressure (mb)
self.ttrace = ma.masked # Parcel trace temperature (C)
self.blayer = ma.masked # Pressure of the bottom of the layer the parcel is lifted (mb)
self.tlayer = ma.masked # Pressure of the top of the layer the parcel is lifted (mb)
self.entrain = 0. # A parcel entrainment setting (not yet implemented)
self.lclpres = ma.masked # Parcel LCL (lifted condensation level) pressure (mb)
self.lclhght = ma.masked # Parcel LCL height (m AGL)
self.lfcpres = ma.masked # Parcel LFC (level of free convection) pressure (mb)
self.lfchght = ma.masked # Parcel LFC height (m AGL)
self.elpres = ma.masked # Parcel EL (equilibrium level) pressure (mb)
self.elhght = ma.masked # Parcel EL height (m AGL)
self.mplpres = ma.masked # Maximum Parcel Level (mb)
self.mplhght = ma.masked # Maximum Parcel Level (m AGL)
self.bplus = ma.masked # Parcel CAPE (J/kg)
self.bminus = ma.masked # Parcel CIN (J/kg)
self.bfzl = ma.masked # Parcel CAPE up to freezing level (J/kg)
self.b3km = ma.masked # Parcel CAPE up to 3 km (J/kg)
self.b6km = ma.masked # Parcel CAPE up to 6 km (J/kg)
self.p0c = ma.masked # Pressure value at 0 C (mb)
self.pm10c = ma.masked # Pressure value at -10 C (mb)
self.pm20c = ma.masked # Pressure value at -20 C (mb)
self.pm30c = ma.masked # Pressure value at -30 C (mb)
self.hght0c = ma.masked # Height value at 0 C (m AGL)
self.hghtm10c = ma.masked # Height value at -10 C (m AGL)
self.hghtm20c = ma.masked # Height value at -20 C (m AGL)
self.hghtm30c = ma.masked # Height value at -30 C (m AGL)
self.wm10c = ma.masked # w velocity at -10 C ?
self.wm20c = ma.masked # w velocity at -20 C ?
self.wm30c = ma.masked # Wet bulb at -30 C ?
self.li5 = ma.masked # Lifted Index at 500 mb (C)
self.li3 = ma.masked # Lifted Index at 300 mb (C)
self.brnshear = ma.masked # Bulk Richardson Number Shear
self.brnu = ma.masked # Bulk Richardson Number U (kts)
self.brnv = ma.masked # Bulk Richardson Number V (kts)
self.brn = ma.masked # Bulk Richardson Number (unitless)
self.limax = ma.masked # Maximum Lifted Index (C)
self.limaxpres = ma.masked # Pressure at Maximum Lifted Index (mb)
self.cap = ma.masked # Cap Strength (C)
self.cappres = ma.masked # Cap strength pressure (mb)
self.bmin = ma.masked # Buoyancy minimum in profile (C)
self.bminpres = ma.masked # Buoyancy minimum pressure (mb)
for kw in kwargs: setattr(self, kw, kwargs.get(kw))
def hgz(prof):
'''
Hail Growth Zone Levels
This function finds the pressure levels for the dendritic
growth zone (from -10 C to -30 C). If either temperature cannot be found,
it is set to be the surface pressure.
Parameters
----------
prof : profile object
Profile Object
Returns
-------
pbot : number
Pressure of the bottom level (mb)
ptop : number
Pressure of the top level (mb)
'''
pbot = temp_lvl(prof, -10)
ptop = temp_lvl(prof, -30)
if not utils.QC(pbot):
pbot = prof.pres[prof.sfc]
if not utils.QC(ptop):
ptop = prof.pres[prof.sfc]
return pbot, ptop
def dgz(prof):
'''
Dendritic Growth Zone Levels
This function finds the pressure levels for the dendritic
growth zone (from -12 C to -17 C). If either temperature cannot be found,
it is set to be the surface pressure.
Parameters
----------
prof : profile object
Profile Object
Returns
-------
pbot : number
Pressure of the bottom level (mb)
ptop : number
Pressure of the top level (mb)
'''
pbot = temp_lvl(prof, -12)
ptop = temp_lvl(prof, -17)
if not utils.QC(pbot):
pbot = prof.pres[prof.sfc]
if not utils.QC(ptop):
ptop = prof.pres[prof.sfc]
return pbot, ptop
def lhp(prof):
'''
Large Hail Parameter
From Johnson and Sugden (2014), EJSSM
.. warning::
This code has not been compared directly against an SPC version.
Parameters
----------
prof : profile object
ConvectiveProfile object
Returns
-------
lhp : number
large hail parameter (unitless)
'''
mag06_shr = utils.KTS2MS(utils.mag(*prof.sfc_6km_shear))
if prof.mupcl.bplus >= 400 and mag06_shr >= 14:
lr75 = prof.lapserate_700_500
zbot, ztop = interp.hght(prof, hgz(prof))
thk_hgz = ztop - zbot
term_a = (((prof.mupcl.bplus - 2000.)/1000.) +\
((3200 - thk_hgz)/500.) +\
((lr75 - 6.5)/2.))
if term_a < 0:
term_a = 0
p_1km, p_3km, p_6km = interp.pres(prof, interp.to_msl(prof, [1000, 3000, 6000]))
shear_el = utils.KTS2MS(utils.mag(*winds.wind_shear(prof, pbot=prof.pres[prof.sfc], ptop=prof.mupcl.elpres)))
grw_el_dir = interp.vec(prof, prof.mupcl.elpres)[0]
grw_36_dir = utils.comp2vec(*winds.mean_wind(prof, pbot=p_3km, ptop=p_6km))[0]
grw_alpha_el = grw_el_dir - grw_36_dir
if grw_alpha_el > 180:
grw_alpha_el = -10
srw_01_dir = utils.comp2vec(*winds.sr_wind(prof, pbot=prof.pres[prof.sfc], ptop=p_1km, stu=prof.srwind[0], stv=prof.srwind[1]))[0]
srw_36_dir = utils.comp2vec(*winds.sr_wind(prof, pbot=p_3km, ptop=p_6km, stu=prof.srwind[0], stv=prof.srwind[1]))[0]
srw_alpha_mid = srw_36_dir - srw_01_dir
term_b = (((shear_el - 25.)/5.) +\
((grw_alpha_el + 5.)/20.) +\
((srw_alpha_mid - 80.)/10.))
if term_b < 0:
term_b = 0
lhp = term_a * term_b + 5
else:
lhp = 0
return lhp
def ship(prof, **kwargs):
'''
Calculate the Sig Hail Parameter (SHIP)
Ryan Jewell (SPC) helped in correcting this equation as the SPC
sounding help page version did not have the correct information
of how SHIP was calculated.
The significant hail parameter (SHIP; SPC 2014) is
an index developed in-house at the SPC. (Johnson and Sugden 2014)
Parameters
----------
prof : profile object
Profile object
mupcl : parcel object, optional
Most Unstable Parcel object
lr75 : float, optional
700 - 500 mb lapse rate (C/km)
h5_temp : float, optional
500 mb temperature (C)
shr06 : float, optional
0-6 km shear (m/s)
frz_lvl : float, optional
freezing level (m)
Returns
-------
ship : number
significant hail parameter (unitless)
'''
mupcl = kwargs.get('mupcl', None)
sfc6shr = kwargs.get('sfc6shr', None)
frz_lvl = kwargs.get('frz_lvl', None)
h5_temp = kwargs.get('h5_temp', None)
lr75 = kwargs.get('lr75', None)
if not mupcl:
try:
mupcl = prof.mupcl
except:
mulplvals = DefineParcel(prof, flag=3, pres=300)
mupcl = cape(prof, lplvals=mulplvals)
mucape = mupcl.bplus
mumr = thermo.mixratio(mupcl.pres, mupcl.dwpc)
if not frz_lvl:
frz_lvl = interp.hght(prof, temp_lvl(prof, 0))
if not h5_temp:
h5_temp = interp.temp(prof, 500.)
if not lr75:
lr75 = lapse_rate(prof, 700., 500., pres=True)
if not sfc6shr:
try:
sfc_6km_shear = prof.sfc_6km_shear
except:
sfc = prof.pres[prof.sfc]
p6km = interp.pres(prof, interp.to_msl(prof, 6000.))
sfc_6km_shear = winds.wind_shear(prof, pbot=sfc, ptop=p6km)
sfc_6km_shear = utils.mag(sfc_6km_shear[0], sfc_6km_shear[1])
shr06 = utils.KTS2MS(sfc_6km_shear)
if shr06 > 27:
shr06 = 27.
elif shr06 < 7:
shr06 = 7.
if mumr > 13.6:
mumr = 13.6
elif mumr < 11.:
mumr = 11.
if h5_temp > -5.5:
h5_temp = -5.5
ship = -1. * (mucape * mumr * lr75 * h5_temp * shr06) / 42000000.
if mucape < 1300:
ship = ship*(mucape/1300.)
if lr75 < 5.8:
ship = ship*(lr75/5.8)
if frz_lvl < 2400:
ship = ship * (frz_lvl/2400.)
return ship
def stp_cin(mlcape, esrh, ebwd, mllcl, mlcinh):
'''
Significant Tornado Parameter (w/CIN)
Formulated using the methodology outlined in [1]_. Used to detect environments where significant tornadoes
are possible within the United States. Uses the effective inflow layer calculations in [3]_ and was created
as an alternative to [2]_.
.. [1] Thompson, R. L., B. T. Smith, J. S. Grams, A. R. Dean, and C. Broyles, 2012: Convective modes for significant severe thunderstorms in the contiguous United States.Part II: Supercell and QLCS tornado environments. Wea. Forecasting, 27, 1136–1154,doi:https://doi.org/10.1175/WAF-D-11-00116.1.
.. [3] Thompson, R. L., C. M. Mead, and R. Edwards, 2007: Effective storm-relative helicity and bulk shear in supercell thunderstorm environments. Wea. Forecasting, 22, 102–115, doi:https://doi.org/10.1175/WAF969.1.
Parameters
----------
mlcape : float
Mixed-layer CAPE from the parcel class (J/kg)
esrh : float
effective storm relative helicity (m2/s2)
ebwd : float
effective bulk wind difference (m/s)
mllcl : float
mixed-layer lifted condensation level (m)
mlcinh : float
mixed-layer convective inhibition (J/kg)
Returns
-------
stp_cin : number
significant tornado parameter (unitless)
See Also
--------
stp_fixed
'''
cape_term = mlcape / 1500.
eshr_term = esrh / 150.
if ebwd < 12.5:
ebwd_term = 0.
elif ebwd > 30.:
ebwd_term = 1.5
else:
ebwd_term = ebwd / 20.
if mllcl < 1000.:
lcl_term = 1.0
elif mllcl > 2000.:
lcl_term = 0.0
else:
lcl_term = ((2000. - mllcl) / 1000.)
if mlcinh > -50:
cinh_term = 1.0
elif mlcinh < -200:
cinh_term = 0
else:
cinh_term = ((mlcinh + 200.) / 150.)
stp_cin = np.maximum(cape_term * eshr_term * ebwd_term * lcl_term * cinh_term, 0)
return stp_cin
def stp_fixed(sbcape, sblcl, srh01, bwd6):
'''
Significant Tornado Parameter (fixed layer)
Formulated using the methodology in [2]_. Used to detect environments where significant tornadoes
are possible within the United States.
.. [2] Thompson, R. L., R. Edwards, J. A. Hart, K. L. Elmore, and P. Markowski, 2003: Close proximity soundings within supercell environments obtained from the Rapid Update Cycle. Wea. Forecasting, 18, 1243–1261, doi:https://doi.org/10.1175/1520-0434(2003)018<1243:CPSWSE>2.0.CO;2
Parameters
----------
sbcape : number
Surface based CAPE from the parcel class (J/kg)
sblcl : number
Surface based lifted condensation level (LCL) (m)
srh01 : number
Surface to 1 km storm relative helicity (m2/s2)
bwd6 : number
Bulk wind difference between 0 to 6 km (m/s)
Returns
-------
stp_fixed : number
signifcant tornado parameter (fixed-layer)
'''
# Calculate SBLCL term
if sblcl < 1000.: # less than 1000. meters
lcl_term = 1.0
elif sblcl > 2000.: # greater than 2000. meters
lcl_term = 0.0
else:
lcl_term = ((2000.-sblcl)/1000.)
# Calculate 6BWD term
if bwd6 > 30.: # greater than 30 m/s
bwd6 = 30
elif bwd6 < 12.5:
bwd6 = 0.0
bwd6_term = bwd6 / 20.
cape_term = sbcape / 1500.
srh_term = srh01 / 150.
stp_fixed = cape_term * lcl_term * srh_term * bwd6_term
return stp_fixed
def scp(mucape, srh, ebwd):
'''
Supercell Composite Parameter
From Thompson et al. 2004, updated from the methodology in [2]_ and uses
the effective inflow layer.
Parameters
----------
prof : profile object
Profile object
mucape : number, optional
Most Unstable CAPE from the parcel class (J/kg) (optional)
srh : number, optional
the effective SRH from the winds.helicity function (m2/s2)
ebwd : number, optional
effective bulk wind difference (m/s)
Returns
-------
scp : number
supercell composite parameter
'''
if ebwd > 20:
ebwd = 20.
elif ebwd < 10:
ebwd = 0.
muCAPE_term = mucape / 1000.
esrh_term = srh / 50.
ebwd_term = ebwd / 20.
scp = muCAPE_term * esrh_term * ebwd_term
return scp
def k_index(prof):
'''
Calculates the K-Index from a profile object
Parameters
----------
prof : profile object
Profile Object
Returns
-------
k_index : number
K-Index
'''
t8 = interp.temp(prof, 850.)
t7 = interp.temp(prof, 700.)
t5 = interp.temp(prof, 500.)
td7 = interp.dwpt(prof, 700.)
td8 = interp.dwpt(prof, 850.)
return t8 - t5 + td8 - (t7 - td7)
def t_totals(prof):
'''
Calculates the Total Totals Index from a profile object
Parameters
----------
prof : profile object
Profile Object
Returns
-------
t_totals : number
Total Totals Index
'''
return c_totals(prof) + v_totals(prof)
def c_totals(prof):
'''
Calculates the Cross Totals Index from a profile object
Parameters
----------
prof : profile object
Profile Object
Returns
-------
c_totals : number
Cross Totals Index
'''
return interp.dwpt(prof, 850.) - interp.temp(prof, 500.)
def v_totals(prof):
'''
Calculates the Vertical Totals Index from a profile object
Parameters
----------
prof : profile object
Profile Object
Returns
-------
v_totals : number
Vertical Totals Index
'''
return interp.temp(prof, 850.) - interp.temp(prof, 500.)
def precip_water(prof, pbot=None, ptop=400, dp=-1, exact=False):
'''
Calculates the precipitable water from a profile object within the
specified layer. The default layer (lower=-1 & upper=-1) is defined to
be surface to 400 hPa.
Parameters
----------
prof : profile object
Profile Object
pbot : number (optional; default surface)
Pressure of the bottom level (hPa)
ptop : number (optional; default 400 hPa)
Pressure of the top level (hPa).
dp : negative integer (optional; default = -1)
The pressure increment for the interpolated sounding
exact : bool (optional; default = False)
Switch to choose between using the exact data (slower) or using
interpolated sounding at 'dp' pressure levels (faster)
Returns
-------
pwat : number,
Precipitable Water (in)
'''
if not pbot: pbot = prof.pres[prof.sfc]
if prof.pres[-1] > ptop:
ptop = prof.pres[-1]
if exact:
ind1 = np.where(pbot > prof.pres)[0].min()
ind2 = np.where(ptop < prof.pres)[0].max()
dwpt1 = interp.dwpt(prof, pbot)
dwpt2 = interp.dwpt(prof, ptop)
mask = ~prof.dwpc.mask[ind1:ind2+1] * ~prof.pres.mask[ind1:ind2+1]
dwpt = np.concatenate([[dwpt1], prof.dwpc[ind1:ind2+1][mask], [dwpt2]])
p = np.concatenate([[pbot], prof.pres[ind1:ind2+1][mask], [ptop]])
else:
dp = -1
p = np.arange(pbot, ptop+dp, dp, dtype=type(pbot))
dwpt = interp.dwpt(prof, p)
w = thermo.mixratio(p, dwpt)
return (((w[:-1]+w[1:])/2 * (p[:-1]-p[1:])) * 0.00040173).sum()
def inferred_temp_adv(prof, dp=-100, lat=35):
'''
Inferred Temperature Advection
SHARP code deduced by Greg Blumberg. Not based on actual SPC code.
Calculates the inferred temperature advection from the surface pressure
and up every 100 mb assuming all winds are geostrophic. The units returned are
in C/hr. If no latitude is specified the function defaults to 35 degrees North.
This code uses Equation 4.1.139 from Bluestein's "Synoptic-Dynamic Meteorology in Midlatitudes (Volume I)"
.. important::
While this code compares well qualitatively to the version at SPC, the SPC output is much larger. Scale
analysis suggests that the values provided by this function are much more reasonable (10 K/day is typical
for synoptic scale values).
Parameters
----------
prof : profile object
Profile object
dp : number, optional
layer size to compute temperature advection over
lat : number, optional
latitude in decimal degrees
Returns
-------
temp_adv : array
an array of temperature advection values (C/hr)
pressure_bounds: array
a 2D array indicating the top and bottom bounds of the temperature advection layers (mb)
'''
if prof.wdir.count() == 0:
return ma.masked, ma.masked
if np.ma.max(prof.pres) <= 100:
return ma.masked, ma.masked
omega = (2. * np.pi) / (86164.)
pres_idx = np.where(prof.pres >= 100.)[0]
pressures = np.arange(prof.pres[prof.get_sfc()], prof.pres[pres_idx][-1], dp, dtype=type(prof.pres[prof.get_sfc()])) # Units: mb
temps = thermo.ctok(interp.temp(prof, pressures))
heights = interp.hght(prof, pressures)
temp_adv = np.empty(len(pressures) - 1)
dirs = interp.vec(prof, pressures)[0]
pressure_bounds = np.empty((len(pressures) - 1, 2))
if utils.QC(lat):
f = 2. * omega * np.sin(np.radians(lat)) # Units: (s**-1)
else:
temp_adv[:] = np.nan
return temp_adv, pressure_bounds
multiplier = (f / 9.81) * (np.pi / 180.) # Units: (s**-1 / (m/s**2)) * (radians/degrees)
for i in range(1, len(pressures)):
bottom_pres = pressures[i-1]
top_pres = pressures[i]
# Get the temperatures from both levels (in Kelvin)
btemp = temps[i-1]
ttemp = temps[i]
# Get the two heights of the top and bottom layer
bhght = heights[i-1] # Units: meters
thght = heights[i] # Units: meters
bottom_wdir = dirs[i-1] # Meteorological degrees (degrees from north)
top_wdir = dirs[i] # same units as top_wdir
# Calculate the average temperature
avg_temp = (ttemp + btemp) * 2.
# Calculate the mean wind between the two levels (this is assumed to be geostrophic)
mean_u, mean_v = winds.mean_wind(prof, pbot=bottom_pres, ptop=top_pres)
mean_wdir, mean_wspd = utils.comp2vec(mean_u, mean_v) # Wind speed is in knots here
mean_wspd = utils.KTS2MS(mean_wspd) # Convert this geostrophic wind speed to m/s
# Here we calculate the change in wind direction with height (thanks to Andrew Mackenzie for help with this)
# The sign of d_theta will dictate whether or not it is warm or cold advection
mod = 180 - bottom_wdir
top_wdir = top_wdir + mod
if top_wdir < 0:
top_wdir = top_wdir + 360
elif top_wdir >= 360:
top_wdir = top_wdir - 360
d_theta = top_wdir - 180.
# Here we calculate t_adv (which is -V_g * del(T) or the local change in temperature term)
# K/s s * rad/m * deg m^2/s^2 K degrees / m
t_adv = multiplier * np.power(mean_wspd,2) * avg_temp * (d_theta / (thght - bhght)) # Units: Kelvin / seconds
# Append the pressure bounds so the person knows the pressure
pressure_bounds[i-1, :] = bottom_pres, top_pres
temp_adv[i-1] = t_adv*60.*60. # Converts Kelvin/seconds to Kelvin/hour (or Celsius/hour)
return temp_adv, pressure_bounds
def temp_lvl(prof, temp, wetbulb=False):
'''
Calculates the level (hPa) of the first occurrence of the specified
temperature.
Parameters
----------
prof : profile object
Profile Object
temp : number
Temperature being searched (C)
wetbulb : boolean
Flag to indicate whether or not the wetbulb profile should be used instead
Returns
-------
First Level of the temperature (hPa) : number
'''
if wetbulb is False:
profile = prof.tmpc
else:
profile = prof.wetbulb
difft = profile - temp
if not np.any(difft <= 0) or not np.any(difft >= 0):
# Temp doesn't occur anywhere; return masked
return ma.masked
elif np.any(difft == 0):
# Temp is one of the data points; don't bother interpolating
return prof.pres[difft == 0][0]
mask = difft.mask | prof.logp.mask
difft = difft[~mask]
profile = profile[~mask]
logp = prof.logp[~mask]
# Find where subsequent values of difft are of opposite sign (i.e. when multiplied together, the result is negative)
ind = np.where((difft[:-1] * difft[1:]) < 0)[0]
try:
ind = ind.min()
except:
ind = ind1[-1]
return np.power(10, np.interp(temp, [profile[ind+1], profile[ind]],
[logp[ind+1], logp[ind]]))