-
Notifications
You must be signed in to change notification settings - Fork 56
/
meteorology.py
1534 lines (1184 loc) · 53.8 KB
/
meteorology.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
import dask.array as da
import numpy as np
import typing
import warnings
import xarray as xr
from .gc_util import _generate_wrapper_docstring
def _dewtemp(
tk: typing.Union[np.ndarray, xr.DataArray, list,
float], rh: typing.Union[np.ndarray, xr.DataArray, list,
float]
) -> typing.Union[np.ndarray, xr.DataArray, list, float]:
"""This function calculates the dew point temperature given temperature and
relative humidity using equations from John Dutton's "Ceaseless Wind" (pp
273-274)
Parameters
----------
tk : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Temperature in Kelvin
rh : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Relative humidity. Must be the same dimensions as ``temperature``
Returns
-------
tdk : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Dewpoint temperature in Kelvin. Same size as input variable temperature
See Also
--------
Related GeoCAT Functions:
:func:`dewtemp`
Related NCL Functions:
`dewtemp_trh <https://www.ncl.ucar.edu/Document/Functions/Built-in/dewtemp_trh.shtml>`__
"""
gc = 461.5 # gas constant for water vapor [j/{kg-k}]
gcx = gc / (1000 * 4.186) # [cal/{g-k}]
lhv = (597.3 - 0.57 * (tk - 273)) / gcx
tdk = tk * lhv / (lhv - tk * np.log(rh * 0.01))
return tdk
def _heat_index(temperature: np.ndarray,
relative_humidity: typing.Union[np.ndarray, xr.DataArray, list,
float],
alternate_coeffs: bool = False) -> np.ndarray:
"""Compute the 'heat index' as calculated by the National Weather Service.
Internal function for heat_index
Parameters
----------
temperature : ndarray
temperature(s) in Fahrenheit
relative_humidity : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
relative humidity as a percentage. Must be the same shape as
``temperature``
alternate_coeffs : bool, optional
flag to use alternate set of coefficients appropriate for
temperatures from 70F to 115F and humidities between 0% and 80%
Returns
-------
heatindex : ndarray
Calculated heat index. Same shape as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`heat_index`,
:func:`_xheat_index`
Related NCL Functions:
`heat_index_nws <https://www.ncl.ucar.edu/Document/Functions/Contributed/heat_index_nws.shtml>`__
"""
# Default coefficients for (t>=80F) and (40<gh<100)
coeffs = [
-42.379, 2.04901523, 10.14333127, -0.22475541, -0.00683783, -0.05481717,
0.00122874, 0.00085282, -0.00000199
]
crit = [80, 40, 100] # [T_low [F], RH_low, RH_high]
# optional flag coefficients for (70F<t<115F) and (0<gh<80)
# within 3F of default coeffs
if alternate_coeffs:
coeffs = [
0.363445176, 0.988622465, 4.777114035, -0.114037667, -0.000850208,
-0.020716198, 0.000687678, 0.000274954, 0.0
]
crit = [70, 0, 80] # [T_low [F], RH_low, RH_high]
# NWS practice
# average Steadman and t
heatindex = (0.5 * (temperature + 61.0 + ((temperature - 68.0) * 1.2) +
(relative_humidity * 0.094)) + temperature) * 0.5
# http://ehp.niehs.nih.gov/1206273/
heatindex = xr.where(temperature < 40, temperature, heatindex)
# if all t values less than critical, return hi
# otherwise perform calculation
eqtype = 0
if not all(temperature.ravel() < crit[0]):
eqtype = 1
heatindex = xr.where(heatindex > crit[0],
_nws_eqn(coeffs, temperature, relative_humidity),
heatindex)
# adjustments
heatindex = xr.where(
np.logical_and(relative_humidity < 13,
np.logical_and(temperature > 80, temperature < 112)),
heatindex - ((13 - relative_humidity) / 4) * np.sqrt(
(17 - abs(temperature - 95)) / 17), heatindex)
heatindex = xr.where(
np.logical_and(relative_humidity > 85,
np.logical_and(temperature > 80, temperature < 87)),
heatindex + ((relative_humidity - 85.0) / 10.0) *
((87.0 - temperature) / 5.0), heatindex)
return heatindex
def _nws_eqn(coeffs, temp, rel_hum):
"""Helper function to compute the heat index.
Internal function for heat_index
Parameters
----------
coeffs : ndarray
coefficients to calculate heat index
temp : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
temperature
rel_hum : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
relative humidity as a percentage. Must be the same shape as
``temperature``
Returns
-------
heatindex : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
Intermediate calculated heat index. Same shape as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`heat_index`,
:func:`_heat_index`
Related NCL Functions:
`heat_index_nws <https://www.ncl.ucar.edu/Document/Functions/Contributed/heat_index_nws.shtml>`__,
"""
heatindex = coeffs[0] \
+ coeffs[1] * temp \
+ coeffs[2] * rel_hum \
+ coeffs[3] * temp * rel_hum \
+ coeffs[4] * temp ** 2 \
+ coeffs[5] * rel_hum ** 2 \
+ coeffs[6] * temp ** 2 * rel_hum \
+ coeffs[7] * temp * rel_hum ** 2 \
+ coeffs[8] * temp ** 2 * rel_hum ** 2
return heatindex
def _relhum(
t: typing.Union[np.ndarray, list, float],
w: typing.Union[np.ndarray, xr.DataArray, list,
float], p: typing.Union[np.ndarray, xr.DataArray, list,
float]) -> np.ndarray:
"""Calculates relative humidity with respect to ice, given temperature,
mixing ratio, and pressure.
"Improved Magnus' Form Approx. of Saturation Vapor pressure"
Oleg A. Alduchov and Robert E. Eskridge
https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml
Parameters
----------
t : ndarray, :obj:`list`, or :obj:`float`
Temperature in Kelvin
w : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
p : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
rh : ndarray
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`relhum_ice`
:func:`relhum_water`
:func:`_xrelhum`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
table = np.asarray([
0.01403, 0.01719, 0.02101, 0.02561, 0.03117, 0.03784, 0.04584, 0.05542,
0.06685, 0.08049, 0.09672, 0.1160, 0.1388, 0.1658, 0.1977, 0.2353,
0.2796, 0.3316, 0.3925, 0.4638, 0.5472, 0.6444, 0.7577, 0.8894, 1.042,
1.220, 1.425, 1.662, 1.936, 2.252, 2.615, 3.032, 3.511, 4.060, 4.688,
5.406, 6.225, 7.159, 8.223, 9.432, 10.80, 12.36, 14.13, 16.12, 18.38,
20.92, 23.80, 27.03, 30.67, 34.76, 39.35, 44.49, 50.26, 56.71, 63.93,
71.98, 80.97, 90.98, 102.1, 114.5, 128.3, 143.6, 160.6, 179.4, 200.2,
223.3, 248.8, 276.9, 307.9, 342.1, 379.8, 421.3, 466.9, 517.0, 572.0,
632.3, 698.5, 770.9, 850.2, 937.0, 1032.0, 1146.6, 1272.0, 1408.1,
1556.7, 1716.9, 1890.3, 2077.6, 2279.6, 2496.7, 2729.8, 2980.0, 3247.8,
3534.1, 3839.8, 4164.8, 4510.5, 4876.9, 5265.1, 5675.2, 6107.8, 6566.2,
7054.7, 7575.3, 8129.4, 8719.2, 9346.50, 10013.0, 10722.0, 11474.0,
12272.0, 13119.0, 14017.0, 14969.0, 15977.0, 17044.0, 18173.0, 19367.0,
20630.0, 21964.0, 23373.0, 24861.0, 26430.0, 28086.0, 29831.0, 31671.0,
33608.0, 35649.0, 37796.0, 40055.0, 42430.0, 44927.0, 47551.0, 50307.0,
53200.0, 56236.0, 59422.0, 62762.0, 66264.0, 69934.0, 73777.0, 77802.0,
82015.0, 86423.0, 91034.0, 95855.0, 100890.0, 106160.0, 111660.0,
117400.0, 123400.0, 129650.0, 136170.0, 142980.0, 150070.0, 157460.0,
165160.0, 173180.0, 181530.0, 190220.0, 199260.0, 208670.0, 218450.0,
228610.0, 239180.0, 250160.0, 261560.0, 273400.0, 285700.0, 298450.0,
311690.0, 325420.0, 339650.0, 354410.0, 369710.0, 385560.0, 401980.0,
418980.0, 436590.0, 454810.0, 473670.0, 493170.0, 513350.0, 534220.0,
555800.0, 578090.0, 601130.0, 624940.0, 649530.0, 674920.0, 701130.0,
728190.0, 756110.0, 784920.0, 814630.0, 845280.0, 876880.0, 909450.0,
943020.0, 977610.0, 1013250.0, 1049940.0, 1087740.0, 1087740.
])
maxtemp = 375.16
mintemp = 173.16
# replace values of t above and below max and min values for temperature
t = np.clip(t, mintemp, maxtemp)
it = (t - mintemp).astype(int)
t2 = mintemp + it
es = (t2 + 1 - t) * table[it] + (t - t2) * table[it + 1]
es = es * 0.1
rh = (w * (p - 0.378 * es) / (0.622 * es)) * 100
# if any value is below 0.0001, set to 0.0001
rh = np.clip(rh, 0.0001, None)
return rh
def _relhum_ice(t: typing.Union[np.ndarray, list, float],
w: typing.Union[np.ndarray, list, float],
p: typing.Union[np.ndarray, list, float]) -> np.ndarray:
"""Calculates relative humidity with respect to ice, given temperature,
mixing ratio, and pressure.
"Improved Magnus' Form Approx. of Saturation Vapor pressure"
Oleg A. Alduchov and Robert E. Eskridge
https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml
Parameters
----------
t : ndarray, :obj:`list`, :obj:`float`
Temperature in Kelvin
w : ndarray, :obj:`list`, :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
p : ndarray, :obj:`list`, :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
rh : ndarray
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`relhum_ice`
:func:`relhum_water`
:func:`_xrelhum`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
# Define data variables
t0 = 273.15
ep = 0.622
onemep = 0.378
es0 = 6.1128
a = 22.571
b = 273.71
est = es0 * np.exp((a * (t - t0)) / ((t - t0) + b))
qst = (ep * est) / ((p * 0.01) - onemep * est)
rh = 100 * (w / qst)
return rh
def _relhum_water(t: typing.Union[np.ndarray, list, float],
w: typing.Union[np.ndarray, list, float],
p: typing.Union[np.ndarray, list, float]) -> np.ndarray:
"""Calculates relative humidity with respect to water, given temperature,
mixing ratio, and pressure.
Definition of mixing ratio if:
- ``es`` - is the saturation mixing ratio
- ``ep`` - is the ratio of the molecular weights of water vapor to dry air
- ``p`` - is the atmospheric pressure
- ``rh`` - is the relative humidity (given as a percent)
.. math::
rh = 100* q / ( (ep*es)/(p-es) )
Parameters
----------
t : ndarray, :obj:`list`, :obj:`float`
Temperature in Kelvin
w : ndarray, :obj:`list`, :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
p : ndarray, :obj:`list`, :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
rh : ndarray
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`relhum_ice`
:func:`relhum_water`
:func:`_xrelhum`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
# Define data variables
t0 = 273.15
ep = 0.622
onemep = 0.378
es0 = 6.1128
a = 17.269
b = 35.86
est = es0 * np.exp((a * (t - t0)) / (t - b))
qst = (ep * est) / ((p * 0.01) - onemep * est)
rh = 100 * (w / qst)
return rh
def _xheat_index(temperature: xr.DataArray,
relative_humidity: xr.DataArray,
alternate_coeffs: bool = False) -> tuple([xr.DataArray, int]):
"""Compute the 'heat index' as calculated by the National Weather Service.
Internal function for heat_index for dask
Parameters
----------
temperature : :class:`xarray.DataArray`
temperature(s) in Fahrenheit
relative_humidity : :class:`xarray.DataArray`
relative humidity as a percentage. Must be the same shape as
``temperature``
alternate_coeffs : bool, optional
flag to use alternate set of coefficients appropriate for
temperatures from 70F to 115F and humidities between 0% and 80%
Returns
-------
heatindex : :class:`xarray.DataArray`
Calculated heat index. Same shape as ``temperature``
eqtype : :class:`int`
version of equations used, for xarray attrs output
See Also
--------
Related GeoCAT Functions:
:func:`heat_index`
:func:`_heat_index`
Related NCL Functions:
`heat_index_nws <https://www.ncl.ucar.edu/Document/Functions/Contributed/heat_index_nws.shtml>`__,
"""
# Default coefficients for (t>=80F) and (40<gh<100)
coeffs = [
-42.379, 2.04901523, 10.14333127, -0.22475541, -0.00683783, -0.05481717,
0.00122874, 0.00085282, -0.00000199
]
crit = [80, 40, 100] # [T_low [F], RH_low, RH_high]
# optional flag coefficients for (70F<t<115F) and (0<gh<80)
# within 3F of default coeffs
if alternate_coeffs:
coeffs = [
0.363445176, 0.988622465, 4.777114035, -0.114037667, -0.000850208,
-0.020716198, 0.000687678, 0.000274954, 0.0
]
crit = [70, 0, 80] # [T_low [F], RH_low, RH_high]
# NWS practice
# average Steadman and t
heatindex = (0.5 * (temperature + 61.0 + ((temperature - 68.0) * 1.2) +
(relative_humidity * 0.094)) + temperature) * 0.5
# http://ehp.niehs.nih.gov/1206273/
heatindex = xr.where(temperature < 40, temperature, heatindex)
# if all t values less than critical, return hi
# otherwise perform calculation
eqtype = 0
if not all(temperature.data.ravel() < crit[0]):
eqtype = 1
heatindex = xr.where(heatindex > crit[0],
_nws_eqn(coeffs, temperature, relative_humidity),
heatindex)
# adjustments
heatindex = xr.where(
np.logical_and(relative_humidity < 13,
np.logical_and(temperature > 80, temperature < 112)),
heatindex - ((13 - relative_humidity) / 4) * np.sqrt(
(17 - abs(temperature - 95)) / 17), heatindex)
heatindex = xr.where(
np.logical_and(relative_humidity > 85,
np.logical_and(temperature > 80, temperature < 87)),
heatindex + ((relative_humidity - 85.0) / 10.0) *
((87.0 - temperature) / 5.0), heatindex)
return heatindex, eqtype
def _xrelhum(t: xr.DataArray, w: xr.DataArray, p: xr.DataArray) -> xr.DataArray:
"""Calculates relative humidity with respect to ice, given temperature,
mixing ratio, and pressure.
"Improved Magnus' Form Approx. of Saturation Vapor pressure"
Oleg A. Alduchov and Robert E. Eskridge
https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml
Parameters
----------
t : :class:`xarray.DataArray`
Temperature in Kelvin
w : :class:`xarray.DataArray`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
p : :class:`xarray.DataArray`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
rh : :class:`xarray.DataArray`
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`relhum_ice`
:func:`relhum_water`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
table = da.from_array([
0.01403, 0.01719, 0.02101, 0.02561, 0.03117, 0.03784, 0.04584, 0.05542,
0.06685, 0.08049, 0.09672, 0.1160, 0.1388, 0.1658, 0.1977, 0.2353,
0.2796, 0.3316, 0.3925, 0.4638, 0.5472, 0.6444, 0.7577, 0.8894, 1.042,
1.220, 1.425, 1.662, 1.936, 2.252, 2.615, 3.032, 3.511, 4.060, 4.688,
5.406, 6.225, 7.159, 8.223, 9.432, 10.80, 12.36, 14.13, 16.12, 18.38,
20.92, 23.80, 27.03, 30.67, 34.76, 39.35, 44.49, 50.26, 56.71, 63.93,
71.98, 80.97, 90.98, 102.1, 114.5, 128.3, 143.6, 160.6, 179.4, 200.2,
223.3, 248.8, 276.9, 307.9, 342.1, 379.8, 421.3, 466.9, 517.0, 572.0,
632.3, 698.5, 770.9, 850.2, 937.0, 1032.0, 1146.6, 1272.0, 1408.1,
1556.7, 1716.9, 1890.3, 2077.6, 2279.6, 2496.7, 2729.8, 2980.0, 3247.8,
3534.1, 3839.8, 4164.8, 4510.5, 4876.9, 5265.1, 5675.2, 6107.8, 6566.2,
7054.7, 7575.3, 8129.4, 8719.2, 9346.50, 10013.0, 10722.0, 11474.0,
12272.0, 13119.0, 14017.0, 14969.0, 15977.0, 17044.0, 18173.0, 19367.0,
20630.0, 21964.0, 23373.0, 24861.0, 26430.0, 28086.0, 29831.0, 31671.0,
33608.0, 35649.0, 37796.0, 40055.0, 42430.0, 44927.0, 47551.0, 50307.0,
53200.0, 56236.0, 59422.0, 62762.0, 66264.0, 69934.0, 73777.0, 77802.0,
82015.0, 86423.0, 91034.0, 95855.0, 100890.0, 106160.0, 111660.0,
117400.0, 123400.0, 129650.0, 136170.0, 142980.0, 150070.0, 157460.0,
165160.0, 173180.0, 181530.0, 190220.0, 199260.0, 208670.0, 218450.0,
228610.0, 239180.0, 250160.0, 261560.0, 273400.0, 285700.0, 298450.0,
311690.0, 325420.0, 339650.0, 354410.0, 369710.0, 385560.0, 401980.0,
418980.0, 436590.0, 454810.0, 473670.0, 493170.0, 513350.0, 534220.0,
555800.0, 578090.0, 601130.0, 624940.0, 649530.0, 674920.0, 701130.0,
728190.0, 756110.0, 784920.0, 814630.0, 845280.0, 876880.0, 909450.0,
943020.0, 977610.0, 1013250.0, 1049940.0, 1087740.0, 1087740.
])
maxtemp = 375.16
mintemp = 173.16
# replace values of t above and below max and min values for temperature
t = t.clip(mintemp, maxtemp)
it = (t - mintemp).astype(int).data
t2 = mintemp + it
it_shape = it.shape
es = (t2 + 1 - t) * table[it.ravel()].reshape(it_shape) + (
t - t2) * table[it.ravel() + 1].reshape(it_shape)
es = es * 0.1
rh = (w * (p - 0.378 * es) / (0.622 * es)) * 100
rh = rh.clip(0.0001, None)
return rh
def dewtemp(
temperature: typing.Union[np.ndarray, xr.DataArray, list, float],
relative_humidity: typing.Union[np.ndarray, xr.DataArray, list, float]
) -> typing.Union[np.ndarray, float]:
"""This function calculates the dew point temperature given temperature and
relative humidity using equations from John Dutton's "Ceaseless Wind" (pp
273-274)
Parameters
----------
temperature : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Temperature in Kelvin
relative_humidity : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Relative humidity. Must be the same dimensions as ``temperature``
Returns
-------
dew_pnt_temp : ndarray or :obj:`float`
Dewpoint temperature in Kelvin. Same size as input variable temperature
See Also
--------
Related GeoCAT Functions:
:func:`_dewtemp`
Related NCL Functions:
`dewtemp_trh <https://www.ncl.ucar.edu/Document/Functions/Built-in/dewtemp_trh.shtml>`__
"""
inputs = [temperature, relative_humidity]
# ensure all inputs same size
if not (np.shape(x) == np.shape(inputs[0]) for x in inputs):
raise ValueError("dewtemp: dimensions of inputs are not the same")
# Get input types
in_types = [type(item) for item in inputs]
if xr.DataArray in in_types:
# check all inputs are xarray.DataArray
if any(x != xr.DataArray for x in in_types):
raise TypeError(
"relhum: if using xarray, all inputs must be xarray")
# call internal computation function
# note: no alternative internal function required for dewtemp
dew_pnt_temp = _dewtemp(temperature, relative_humidity)
# set xarray attributes
dew_pnt_temp.attrs['long_name'] = 'dew point temperature'
dew_pnt_temp.attrs['units'] = 'Kelvin'
else:
# ensure in numpy array for function call
temperature = np.asarray(temperature)
relative_humidity = np.asarray(relative_humidity)
# function call for non-dask/xarray
dew_pnt_temp = _dewtemp(temperature, relative_humidity)
return dew_pnt_temp
def heat_index(
temperature: typing.Union[np.ndarray, xr.DataArray, list, float],
relative_humidity: typing.Union[np.ndarray, xr.DataArray, list, float],
alternate_coeffs: bool = False
) -> typing.Union[np.ndarray, xr.DataArray]:
"""Compute the 'heat index' as calculated by the National Weather Service.
The heat index calculation in this function is described at:
https://www.wpc.ncep.noaa.gov/html/heatindex_equation.shtml
The 'Heat Index' is a measure of how hot weather "feels" to the body. The combination of temperature and humidity
produce an "apparent temperature" or the temperature the body "feels". The returned values are for shady
locations only. Exposure to full sunshine can increase heat index values by up to 15°F. Also, strong winds,
particularly with very hot, dry air, can be extremely hazardous as the wind adds heat to the body
The computation of the heat index is a refinement of a result obtained by multiple regression analysis carried
out by Lans P. Rothfusz and described in a 1990 National Weather Service (NWS) Technical Attachment (SR 90-23).
All values less that 40F/4.4C/277.65K are set to the ambient temperature.
In practice, the Steadman formula is computed first and the result averaged with the temperature. If this heat
index value is 80 degrees F or higher, the full regression equation along with any adjustment as described above
is applied. If the ambient temperature is less the 40F (4.4C/277.65K), the heat index is set to the ambient
temperature.
Parameters
----------
temperature : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
temperature(s) in Fahrenheit
relative_humidity : ndarray, :class:`xarray.DataArray`, :class:`list`, :class:`float`
relative humidity as a percentage. Must be the same shape as
``temperature``
alternate_coeffs : bool, optional
flag to use alternate set of coefficients appropriate for
temperatures from 70F to 115F and humidities between 0% and 80%
Returns
-------
heatindex : ndarray, :class:`xarray.DataArray`
Calculated heat index. Same shape as ``temperature``
Examples
--------
>>> import numpy as np
>>> import geocat.comp
>>> t = np.array([104, 100, 92])
>>> rh = np.array([55, 65, 60])
>>> hi = heat_index(t,rh)
>>> hi
array([137.36135724, 135.8679973 , 104.68441864])
See Also
--------
Related GeoCAT Functions:
:func:`_heat_index`
:func:`_xheat_index`
Related NCL Functions:
`heat_index_nws <https://www.ncl.ucar.edu/Document/Functions/Contributed/heat_index_nws.shtml>`__
"""
inputs = [temperature, relative_humidity]
# ensure all inputs same size
if not (np.shape(x) == np.shape(inputs[0]) for x in inputs):
raise ValueError("heat_index: dimensions of inputs are not the same")
# Get input types
in_types = [type(item) for item in inputs]
# run dask compatible version if input is xarray
if xr.DataArray in in_types:
# check all inputs are xarray.DataArray
if not all(x == xr.DataArray for x in in_types):
raise TypeError(
"heat_index: if using xarray, all inputs must be xarray")
# input validation on relative humidity
if any(relative_humidity.data.ravel() < 0) or any(
relative_humidity.data.ravel() > 100):
raise ValueError('heat_index: invalid values for relative humidity')
# Check if relative humidity fractional
if all(relative_humidity.data.ravel() < 1):
warnings.warn(
"WARNING: rh must be %, not fractional; All rh are < 1")
# call internal computation function
heatindex, eqtype = _xheat_index(temperature, relative_humidity,
alternate_coeffs)
# set xarray attributes
heatindex.attrs['long_name'] = "heat index: NWS"
heatindex.attrs['units'] = "F"
heatindex.attrs[
'www'] = "https://www.wpc.ncep.noaa.gov/html/heatindex_equation.shtml"
heatindex.attrs['info'] = "appropriate for shady locations with no wind"
if eqtype == 1:
heatindex.attrs[
'tag'] = "NCL: heat_index_nws; (Steadman+t)*0.5 and Rothfusz"
else:
heatindex.attrs['tag'] = "NCL: heat_index_nws; (Steadman+t)*0.5"
else:
# ensure in numpy array for function call
temperature = np.atleast_1d(temperature)
relative_humidity = np.atleast_1d(relative_humidity)
# input validation on relative humidity
if any(relative_humidity.ravel() < 0) or any(
relative_humidity.ravel() > 100):
raise ValueError('heat_index: invalid values for relative humidity')
# Check if relative humidity fractional
if all(relative_humidity.ravel() < 1):
warnings.warn(
"WARNING: rh must be %, not fractional; All rh are < 1")
# function call for non-dask/xarray
heatindex = _heat_index(temperature, relative_humidity,
alternate_coeffs)
return heatindex
def relhum(
temperature: typing.Union[np.ndarray, xr.DataArray, list, float],
mixing_ratio: typing.Union[np.ndarray, xr.DataArray, list, float],
pressure: typing.Union[np.ndarray, xr.DataArray, list, float]
) -> typing.Union[np.ndarray, xr.DataArray]:
"""This function calculates the relative humidity given temperature, mixing
ratio, and pressure.
"Improved Magnus' Form Approx. of Saturation Vapor pressure"
Oleg A. Alduchov and Robert E. Eskridge
https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml
Parameters
----------
temperature : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Temperature in Kelvin
mixing_ratio : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
pressure : ndarray, :class:`xarray.DataArray`, :obj:`list`, or :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
relative_humidity : ndarray or :class:`xarray.DataArray`
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`_xrelhum`
:func:`relhum_ice`
:func:`relhum_water`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
inputs = [temperature, mixing_ratio, pressure]
# ensure all inputs same size
if not (np.shape(x) == np.shape(inputs[0]) for x in inputs):
raise ValueError("relhum: dimensions of inputs are not the same")
# Get input types
in_types = [type(item) for item in inputs]
# run dask compatible version if input is xarray
if xr.DataArray in in_types:
# check all inputs are xarray.DataArray
if not all(x == xr.DataArray for x in in_types):
raise TypeError(
"relhum: if using xarray, all inputs must be xarray")
# call internal computation function
relative_humidity = _xrelhum(temperature, mixing_ratio, pressure)
# set xarray attributes
relative_humidity.attrs['long_name'] = "relative humidity"
relative_humidity.attrs['units'] = 'percentage'
relative_humidity.attrs[
'info'] = 'https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml'
else:
# ensure in numpy array for function call
temperature = np.asarray(temperature)
mixing_ratio = np.asarray(mixing_ratio)
pressure = np.asarray(pressure)
# function call for non-dask/xarray
relative_humidity = _relhum(temperature, mixing_ratio, pressure)
return relative_humidity
def relhum_ice(temperature: typing.Union[np.ndarray, list, float],
mixing_ratio: typing.Union[np.ndarray, list, float],
pressure: typing.Union[np.ndarray, list, float]) -> np.ndarray:
"""Calculates relative humidity with respect to ice, given temperature,
mixing ratio, and pressure.
"Improved Magnus' Form Approx. of Saturation Vapor pressure"
Oleg A. Alduchov and Robert E. Eskridge
https://journals.ametsoc.org/view/journals/apme/35/4/1520-0450_1996_035_0601_imfaos_2_0_co_2.xml
Parameters
----------
temperature : ndarray, :obj:`list`, or :obj:`float`
Temperature in Kelvin
mixing_ratio : ndarray, :obj:`list`, or :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
pressure : ndarray, :obj:`list`, or :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
relative_humidity : ndarray
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`_xrelhum`
:func:`relhum_water`
:func:`_relhum_ice`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
# If xarray input, pull data and store metadata
x_out = False
if isinstance(temperature, xr.DataArray):
x_out = True
save_dims = temperature.dims
save_coords = temperature.coords
save_attrs = temperature.attrs
# ensure in numpy array for function call
temperature = np.asarray(temperature)
mixing_ratio = np.asarray(mixing_ratio)
pressure = np.asarray(pressure)
# ensure all inputs same size
if np.shape(temperature) != np.shape(mixing_ratio) or np.shape(
temperature) != np.shape(pressure):
raise ValueError("relhum_ice: dimensions of inputs are not the same")
relative_humidity = _relhum_ice(temperature, mixing_ratio, pressure)
# output as xarray if input as xarray
if x_out:
relative_humidity = xr.DataArray(data=relative_humidity,
coords=save_coords,
dims=save_dims,
attrs=save_attrs)
return relative_humidity
def relhum_water(temperature: typing.Union[np.ndarray, list, float],
mixing_ratio: typing.Union[np.ndarray, list, float],
pressure: typing.Union[np.ndarray, list, float]) -> np.ndarray:
"""Calculates relative humidity with respect to water, given temperature,
mixing ratio, and pressure.
Definition of mixing ratio if:
- `es` - is the saturation mixing ratio
- `ep` - is the ratio of the molecular weights of water vapor to dry air
- `p` - is the atmospheric pressure
- `rh` - is the relative humidity (given as a percent)
.. math::
rh = 100 q / ( (ep*es)/(p-es) )
Parameters
----------
temperature : ndarray, :obj:`list`, or :obj:`float`
Temperature in Kelvin
mixing_ratio : ndarray, :obj:`list`, or :obj:`float`
Mixing ratio in kg/kg. Must have the same dimensions as ``temperature``
pressure : ndarray, :obj:`list`, or :obj:`float`
Pressure in Pa. Must have the same dimensions as ``temperature``
Returns
-------
relative_humidity : ndarray
Relative humidity. Will have the same dimensions as ``temperature``
See Also
--------
Related GeoCAT Functions:
:func:`relhum`
:func:`_xrelhum`
:func:`relhum_ice`
:func:`relhum_water`
Related NCL Functions:
`relhum <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum.shtml>`__,
`relhum_ttd <https://www.ncl.ucar.edu/Document/Functions/Contributed/relhum_ttd.shtml>`__,
`relhum_ice <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_ice.shtml>`__,
`relhum_water <https://www.ncl.ucar.edu/Document/Functions/Built-in/relhum_water.shtml>`__
"""
# If xarray input, pull data and store metadata
x_out = False
if isinstance(temperature, xr.DataArray):
x_out = True
save_dims = temperature.dims
save_coords = temperature.coords
save_attrs = temperature.attrs
# ensure in numpy array for function call
temperature = np.asarray(temperature)
mixing_ratio = np.asarray(mixing_ratio)
pressure = np.asarray(pressure)
# ensure all inputs same size
if np.shape(temperature) != np.shape(mixing_ratio) or np.shape(
temperature) != np.shape(pressure):
raise ValueError("relhum_water: dimensions of inputs are not the same")
relative_humidity = _relhum_water(temperature, mixing_ratio, pressure)
# output as xarray if input as xarray
if x_out:
relative_humidity = xr.DataArray(data=relative_humidity,
coords=save_coords,
dims=save_dims,
attrs=save_attrs)
return relative_humidity
def max_daylight(
jday: typing.Union[np.ndarray, xr.DataArray, list,
float], lat: typing.Union[np.ndarray, xr.DataArray, list,
float]
) -> typing.Union[np.ndarray, xr.DataArray, float]: