-
Notifications
You must be signed in to change notification settings - Fork 100
/
climate.py
1448 lines (1206 loc) · 53.4 KB
/
climate.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
"""Climate data and mass-balance computations"""
# Built ins
import logging
import os
import datetime
import warnings
# External libs
import numpy as np
import netCDF4
import pandas as pd
import xarray as xr
from scipy import stats
from scipy import optimize as optimization
# Optional libs
try:
import salem
except ImportError:
pass
# Locals
from oggm import cfg
from oggm import utils
from oggm.core import centerlines
from oggm import entity_task, global_task
from oggm.exceptions import MassBalanceCalibrationError, InvalidParamsError
# Module logger
log = logging.getLogger(__name__)
@entity_task(log, writes=['climate_monthly', 'climate_info'])
def process_custom_climate_data(gdir):
"""Processes and writes the climate data from a user-defined climate file.
The input file must have a specific format (see
https://github.com/OGGM/oggm-sample-data test-files/histalp_merged_hef.nc
for an example).
This is the way OGGM used to do it for HISTALP before it got automatised.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
the glacier directory to process
"""
if not (('climate_file' in cfg.PATHS) and
os.path.exists(cfg.PATHS['climate_file'])):
raise InvalidParamsError('Custom climate file not found')
if cfg.PARAMS['baseline_climate'] not in ['', 'CUSTOM']:
raise InvalidParamsError("When using custom climate data please set "
"PARAMS['baseline_climate'] to an empty "
"string or `CUSTOM`. Note also that you can "
"now use the `process_histalp_data` task for "
"automated HISTALP data processing.")
# read the file
fpath = cfg.PATHS['climate_file']
nc_ts = salem.GeoNetcdf(fpath)
# set temporal subset for the ts data (hydro years)
sm = cfg.PARAMS['hydro_month_' + gdir.hemisphere]
em = sm - 1 if (sm > 1) else 12
yrs = nc_ts.time.year
y0, y1 = yrs[0], yrs[-1]
if cfg.PARAMS['baseline_y0'] != 0:
y0 = cfg.PARAMS['baseline_y0']
if cfg.PARAMS['baseline_y1'] != 0:
y1 = cfg.PARAMS['baseline_y1']
nc_ts.set_period(t0='{}-{:02d}-01'.format(y0, sm),
t1='{}-{:02d}-01'.format(y1, em))
time = nc_ts.time
ny, r = divmod(len(time), 12)
if r != 0:
raise InvalidParamsError('Climate data should be full years')
# Units
assert nc_ts._nc.variables['hgt'].units.lower() in ['m', 'meters', 'meter',
'metres', 'metre']
assert nc_ts._nc.variables['temp'].units.lower() in ['degc', 'degrees',
'degree', 'c']
assert nc_ts._nc.variables['prcp'].units.lower() in ['kg m-2', 'l m-2',
'mm', 'millimeters',
'millimeter']
# geoloc
lon = nc_ts._nc.variables['lon'][:]
lat = nc_ts._nc.variables['lat'][:]
ilon = np.argmin(np.abs(lon - gdir.cenlon))
ilat = np.argmin(np.abs(lat - gdir.cenlat))
ref_pix_lon = lon[ilon]
ref_pix_lat = lat[ilat]
# read the data
temp = nc_ts.get_vardata('temp')
prcp = nc_ts.get_vardata('prcp')
hgt = nc_ts.get_vardata('hgt')
ttemp = temp[:, ilat-1:ilat+2, ilon-1:ilon+2]
itemp = ttemp[:, 1, 1]
thgt = hgt[ilat-1:ilat+2, ilon-1:ilon+2]
ihgt = thgt[1, 1]
thgt = thgt.flatten()
iprcp = prcp[:, ilat, ilon]
nc_ts.close()
# Should we compute the gradient?
use_grad = cfg.PARAMS['temp_use_local_gradient']
igrad = None
if use_grad:
igrad = np.zeros(len(time)) * np.NaN
for t, loct in enumerate(ttemp):
slope, _, _, p_val, _ = stats.linregress(thgt,
loct.flatten())
igrad[t] = slope if (p_val < 0.01) else np.NaN
gdir.write_monthly_climate_file(time, iprcp, itemp, ihgt,
ref_pix_lon, ref_pix_lat,
gradient=igrad)
# metadata
out = {'baseline_climate_source': fpath,
'baseline_hydro_yr_0': y0+1,
'baseline_hydro_yr_1': y1}
gdir.write_pickle(out, 'climate_info')
@entity_task(log, writes=['climate_monthly', 'climate_info'])
def process_cru_data(gdir):
"""Processes and writes the CRU baseline climate data for this glacier.
Interpolates the CRU TS data to the high-resolution CL2 climatologies
(provided with OGGM) and writes everything to a NetCDF file.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
the glacier directory to process
"""
if cfg.PATHS.get('climate_file', None):
warnings.warn("You seem to have set a custom climate file for this "
"run, but are using the default CRU climate "
"file instead.")
if cfg.PARAMS['baseline_climate'] != 'CRU':
raise InvalidParamsError("cfg.PARAMS['baseline_climate'] should be "
"set to CRU")
# read the climatology
clfile = utils.get_cru_cl_file()
ncclim = salem.GeoNetcdf(clfile)
# and the TS data
nc_ts_tmp = salem.GeoNetcdf(utils.get_cru_file('tmp'), monthbegin=True)
nc_ts_pre = salem.GeoNetcdf(utils.get_cru_file('pre'), monthbegin=True)
# set temporal subset for the ts data (hydro years)
sm = cfg.PARAMS['hydro_month_' + gdir.hemisphere]
em = sm - 1 if (sm > 1) else 12
yrs = nc_ts_pre.time.year
y0, y1 = yrs[0], yrs[-1]
if cfg.PARAMS['baseline_y0'] != 0:
y0 = cfg.PARAMS['baseline_y0']
if cfg.PARAMS['baseline_y1'] != 0:
y1 = cfg.PARAMS['baseline_y1']
nc_ts_tmp.set_period(t0='{}-{:02d}-01'.format(y0, sm),
t1='{}-{:02d}-01'.format(y1, em))
nc_ts_pre.set_period(t0='{}-{:02d}-01'.format(y0, sm),
t1='{}-{:02d}-01'.format(y1, em))
time = nc_ts_pre.time
ny, r = divmod(len(time), 12)
assert r == 0
lon = gdir.cenlon
lat = gdir.cenlat
# This is guaranteed to work because I prepared the file (I hope)
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
# get climatology data
loc_hgt = ncclim.get_vardata('elev')
loc_tmp = ncclim.get_vardata('temp')
loc_pre = ncclim.get_vardata('prcp')
loc_lon = ncclim.get_vardata('lon')
loc_lat = ncclim.get_vardata('lat')
# see if the center is ok
if not np.isfinite(loc_hgt[1, 1]):
# take another candidate where finite
isok = np.isfinite(loc_hgt)
# wait: some areas are entirely NaNs, make the subset larger
_margin = 1
while not np.any(isok):
_margin += 1
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=_margin)
loc_hgt = ncclim.get_vardata('elev')
isok = np.isfinite(loc_hgt)
if _margin > 1:
log.debug('(%s) I had to look up for far climate pixels: %s',
gdir.rgi_id, _margin)
# Take the first candidate (doesn't matter which)
lon, lat = ncclim.grid.ll_coordinates
lon = lon[isok][0]
lat = lat[isok][0]
# Resubset
ncclim.set_subset()
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
loc_hgt = ncclim.get_vardata('elev')
loc_tmp = ncclim.get_vardata('temp')
loc_pre = ncclim.get_vardata('prcp')
loc_lon = ncclim.get_vardata('lon')
loc_lat = ncclim.get_vardata('lat')
assert np.isfinite(loc_hgt[1, 1])
isok = np.isfinite(loc_hgt)
hgt_f = loc_hgt[isok].flatten()
assert len(hgt_f) > 0.
# Should we compute the gradient?
use_grad = cfg.PARAMS['temp_use_local_gradient']
ts_grad = None
if use_grad and len(hgt_f) >= 5:
ts_grad = np.zeros(12) * np.NaN
for i in range(12):
loc_tmp_mth = loc_tmp[i, ...][isok].flatten()
slope, _, _, p_val, _ = stats.linregress(hgt_f, loc_tmp_mth)
ts_grad[i] = slope if (p_val < 0.01) else np.NaN
# convert to a timeseries and hydrological years
ts_grad = ts_grad.tolist()
ts_grad = ts_grad[em:] + ts_grad[0:em]
ts_grad = np.asarray(ts_grad * ny)
# maybe this will throw out of bounds warnings
nc_ts_tmp.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
nc_ts_pre.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
# compute monthly anomalies
# of temp
ts_tmp = nc_ts_tmp.get_vardata('tmp', as_xarray=True)
ts_tmp_avg = ts_tmp.sel(time=slice('1961-01-01', '1990-12-01'))
ts_tmp_avg = ts_tmp_avg.groupby('time.month').mean(dim='time')
ts_tmp = ts_tmp.groupby('time.month') - ts_tmp_avg
# of precip
ts_pre = nc_ts_pre.get_vardata('pre', as_xarray=True)
ts_pre_avg = ts_pre.sel(time=slice('1961-01-01', '1990-12-01'))
ts_pre_avg = ts_pre_avg.groupby('time.month').mean(dim='time')
ts_pre_ano = ts_pre.groupby('time.month') - ts_pre_avg
# scaled anomalies is the default. Standard anomalies above
# are used later for where ts_pre_avg == 0
ts_pre = ts_pre.groupby('time.month') / ts_pre_avg
# interpolate to HR grid
if np.any(~np.isfinite(ts_tmp[:, 1, 1])):
# Extreme case, middle pix is not valid
# take any valid pix from the 3*3 (and hope there's one)
found_it = False
for idi in range(2):
for idj in range(2):
if np.all(np.isfinite(ts_tmp[:, idj, idi])):
ts_tmp[:, 1, 1] = ts_tmp[:, idj, idi]
ts_pre[:, 1, 1] = ts_pre[:, idj, idi]
ts_pre_ano[:, 1, 1] = ts_pre_ano[:, idj, idi]
found_it = True
if not found_it:
msg = '({}) there is no climate data'.format(gdir.rgi_id)
raise MassBalanceCalibrationError(msg)
elif np.any(~np.isfinite(ts_tmp)):
# maybe the side is nan, but we can do nearest
ts_tmp = ncclim.grid.map_gridded_data(ts_tmp.values, nc_ts_tmp.grid,
interp='nearest')
ts_pre = ncclim.grid.map_gridded_data(ts_pre.values, nc_ts_pre.grid,
interp='nearest')
ts_pre_ano = ncclim.grid.map_gridded_data(ts_pre_ano.values,
nc_ts_pre.grid,
interp='nearest')
else:
# We can do bilinear
ts_tmp = ncclim.grid.map_gridded_data(ts_tmp.values, nc_ts_tmp.grid,
interp='linear')
ts_pre = ncclim.grid.map_gridded_data(ts_pre.values, nc_ts_pre.grid,
interp='linear')
ts_pre_ano = ncclim.grid.map_gridded_data(ts_pre_ano.values,
nc_ts_pre.grid,
interp='linear')
# take the center pixel and add it to the CRU CL clim
# for temp
loc_tmp = xr.DataArray(loc_tmp[:, 1, 1], dims=['month'],
coords={'month': ts_tmp_avg.month})
ts_tmp = xr.DataArray(ts_tmp[:, 1, 1], dims=['time'],
coords={'time': time})
ts_tmp = ts_tmp.groupby('time.month') + loc_tmp
# for prcp
loc_pre = xr.DataArray(loc_pre[:, 1, 1], dims=['month'],
coords={'month': ts_pre_avg.month})
ts_pre = xr.DataArray(ts_pre[:, 1, 1], dims=['time'],
coords={'time': time})
ts_pre_ano = xr.DataArray(ts_pre_ano[:, 1, 1], dims=['time'],
coords={'time': time})
# scaled anomalies
ts_pre = ts_pre.groupby('time.month') * loc_pre
# standard anomalies
ts_pre_ano = ts_pre_ano.groupby('time.month') + loc_pre
# Correct infinite values with standard anomalies
ts_pre.values = np.where(np.isfinite(ts_pre.values),
ts_pre.values,
ts_pre_ano.values)
# The last step might create negative values (unlikely). Clip them
ts_pre.values = utils.clip_min(ts_pre.values, 0)
# done
loc_hgt = loc_hgt[1, 1]
loc_lon = loc_lon[1]
loc_lat = loc_lat[1]
assert np.isfinite(loc_hgt)
assert np.all(np.isfinite(ts_pre.values))
assert np.all(np.isfinite(ts_tmp.values))
gdir.write_monthly_climate_file(time, ts_pre.values, ts_tmp.values,
loc_hgt, loc_lon, loc_lat,
gradient=ts_grad)
source = nc_ts_tmp._nc.title[:10]
ncclim._nc.close()
nc_ts_tmp._nc.close()
nc_ts_pre._nc.close()
# metadata
out = {'baseline_climate_source': source,
'baseline_hydro_yr_0': y0+1,
'baseline_hydro_yr_1': y1}
gdir.write_pickle(out, 'climate_info')
@entity_task(log, writes=['climate_monthly', 'climate_info'])
def process_dummy_cru_file(gdir, sigma_temp=2, sigma_prcp=0.5, seed=None):
"""Create a simple baseline climate file for this glacier - for testing!
This simply reproduces the climatology with a little randomness in it.
TODO: extend the functionality by allowing a monthly varying sigma
Parameters
----------
gdir : GlacierDirectory
the glacier directory
sigma_temp : float
the standard deviation of the random timeseries (set to 0 for constant
ts)
sigma_prcp : float
the standard deviation of the random timeseries (set to 0 for constant
ts)
seed : int
the RandomState seed
"""
# read the climatology
clfile = utils.get_cru_cl_file()
ncclim = salem.GeoNetcdf(clfile)
# set temporal subset for the ts data (hydro years)
sm = cfg.PARAMS['hydro_month_' + gdir.hemisphere]
em = sm - 1 if (sm > 1) else 12
y0, y1 = 1901, 2018
if cfg.PARAMS['baseline_y0'] != 0:
y0 = cfg.PARAMS['baseline_y0']
if cfg.PARAMS['baseline_y1'] != 0:
y1 = cfg.PARAMS['baseline_y1']
time = pd.date_range(start='{}-{:02d}-01'.format(y0, sm),
end='{}-{:02d}-01'.format(y1, em),
freq='MS')
ny, r = divmod(len(time), 12)
assert r == 0
lon = gdir.cenlon
lat = gdir.cenlat
# This is guaranteed to work because I prepared the file (I hope)
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
# get climatology data
loc_hgt = ncclim.get_vardata('elev')
loc_tmp = ncclim.get_vardata('temp')
loc_pre = ncclim.get_vardata('prcp')
loc_lon = ncclim.get_vardata('lon')
loc_lat = ncclim.get_vardata('lat')
# see if the center is ok
if not np.isfinite(loc_hgt[1, 1]):
# take another candidate where finite
isok = np.isfinite(loc_hgt)
# wait: some areas are entirely NaNs, make the subset larger
_margin = 1
while not np.any(isok):
_margin += 1
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=_margin)
loc_hgt = ncclim.get_vardata('elev')
isok = np.isfinite(loc_hgt)
if _margin > 1:
log.debug('(%s) I had to look up for far climate pixels: %s',
gdir.rgi_id, _margin)
# Take the first candidate (doesn't matter which)
lon, lat = ncclim.grid.ll_coordinates
lon = lon[isok][0]
lat = lat[isok][0]
# Resubset
ncclim.set_subset()
ncclim.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
loc_hgt = ncclim.get_vardata('elev')
loc_tmp = ncclim.get_vardata('temp')
loc_pre = ncclim.get_vardata('prcp')
loc_lon = ncclim.get_vardata('lon')
loc_lat = ncclim.get_vardata('lat')
assert np.isfinite(loc_hgt[1, 1])
isok = np.isfinite(loc_hgt)
hgt_f = loc_hgt[isok].flatten()
assert len(hgt_f) > 0.
# Should we compute the gradient?
use_grad = cfg.PARAMS['temp_use_local_gradient']
ts_grad = None
if use_grad and len(hgt_f) >= 5:
ts_grad = np.zeros(12) * np.NaN
for i in range(12):
loc_tmp_mth = loc_tmp[i, ...][isok].flatten()
slope, _, _, p_val, _ = stats.linregress(hgt_f, loc_tmp_mth)
ts_grad[i] = slope if (p_val < 0.01) else np.NaN
# convert to a timeseries and hydrological years
ts_grad = ts_grad.tolist()
ts_grad = ts_grad[em:] + ts_grad[0:em]
ts_grad = np.asarray(ts_grad * ny)
# Make DataArrays
rng = np.random.RandomState(seed)
loc_tmp = xr.DataArray(loc_tmp[:, 1, 1], dims=['month'],
coords={'month': np.arange(1, 13)})
ts_tmp = rng.randn(len(time)) * sigma_temp
ts_tmp = xr.DataArray(ts_tmp, dims=['time'],
coords={'time': time})
loc_pre = xr.DataArray(loc_pre[:, 1, 1], dims=['month'],
coords={'month': np.arange(1, 13)})
ts_pre = utils.clip_min(rng.randn(len(time)) * sigma_prcp + 1, 0)
ts_pre = xr.DataArray(ts_pre, dims=['time'],
coords={'time': time})
# Create the time series
ts_tmp = ts_tmp.groupby('time.month') + loc_tmp
ts_pre = ts_pre.groupby('time.month') * loc_pre
# done
loc_hgt = loc_hgt[1, 1]
loc_lon = loc_lon[1]
loc_lat = loc_lat[1]
assert np.isfinite(loc_hgt)
gdir.write_monthly_climate_file(time, ts_pre.values, ts_tmp.values,
loc_hgt, loc_lon, loc_lat,
gradient=ts_grad)
source = 'CRU CL2 and some randomness'
ncclim._nc.close()
# metadata
out = {'baseline_climate_source': source,
'baseline_hydro_yr_0': y0+1,
'baseline_hydro_yr_1': y1}
gdir.write_pickle(out, 'climate_info')
@entity_task(log, writes=['climate_monthly', 'climate_info'])
def process_histalp_data(gdir):
"""Processes and writes the HISTALP baseline climate data for this glacier.
Extracts the nearest timeseries and writes everything to a NetCDF file.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
the glacier directory to process
"""
if cfg.PATHS.get('climate_file', None):
warnings.warn("You seem to have set a custom climate file for this "
"run, but are using the default HISTALP climate file "
"instead.")
if cfg.PARAMS['baseline_climate'] != 'HISTALP':
raise InvalidParamsError("cfg.PARAMS['baseline_climate'] should be "
"set to HISTALP.")
# read the time out of the pure netcdf file
ft = utils.get_histalp_file('tmp')
fp = utils.get_histalp_file('pre')
with utils.ncDataset(ft) as nc:
vt = nc.variables['time']
assert vt[0] == 0
assert vt[-1] == vt.shape[0] - 1
t0 = vt.units.split(' since ')[1][:7]
time_t = pd.date_range(start=t0, periods=vt.shape[0], freq='MS')
with utils.ncDataset(fp) as nc:
vt = nc.variables['time']
assert vt[0] == 0.5
assert vt[-1] == vt.shape[0] - .5
t0 = vt.units.split(' since ')[1][:7]
time_p = pd.date_range(start=t0, periods=vt.shape[0], freq='MS')
# Now open with salem
nc_ts_tmp = salem.GeoNetcdf(ft, time=time_t)
nc_ts_pre = salem.GeoNetcdf(fp, time=time_p)
# set temporal subset for the ts data (hydro years)
# the reference time is given by precip, which is shorter
sm = cfg.PARAMS['hydro_month_nh']
em = sm - 1 if (sm > 1) else 12
yrs = nc_ts_pre.time.year
y0, y1 = yrs[0], yrs[-1]
if cfg.PARAMS['baseline_y0'] != 0:
y0 = cfg.PARAMS['baseline_y0']
if cfg.PARAMS['baseline_y1'] != 0:
y1 = cfg.PARAMS['baseline_y1']
nc_ts_tmp.set_period(t0='{}-{:02d}-01'.format(y0, sm),
t1='{}-{:02d}-01'.format(y1, em))
nc_ts_pre.set_period(t0='{}-{:02d}-01'.format(y0, sm),
t1='{}-{:02d}-01'.format(y1, em))
time = nc_ts_pre.time
ny, r = divmod(len(time), 12)
assert r == 0
# Units
assert nc_ts_tmp._nc.variables['HSURF'].units.lower() in ['m', 'meters',
'meter',
'metres',
'metre']
assert nc_ts_tmp._nc.variables['T_2M'].units.lower() in ['degc', 'degrees',
'degrees celcius',
'degree', 'c']
assert nc_ts_pre._nc.variables['TOT_PREC'].units.lower() in ['kg m-2',
'l m-2', 'mm',
'millimeters',
'millimeter']
# geoloc
lon = gdir.cenlon
lat = gdir.cenlat
nc_ts_tmp.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
nc_ts_pre.set_subset(corners=((lon, lat), (lon, lat)), margin=1)
# read the data
temp = nc_ts_tmp.get_vardata('T_2M')
prcp = nc_ts_pre.get_vardata('TOT_PREC')
hgt = nc_ts_tmp.get_vardata('HSURF')
ref_lon = nc_ts_tmp.get_vardata('lon')
ref_lat = nc_ts_tmp.get_vardata('lat')
source = nc_ts_tmp._nc.title[:7]
nc_ts_tmp._nc.close()
nc_ts_pre._nc.close()
# Should we compute the gradient?
use_grad = cfg.PARAMS['temp_use_local_gradient']
igrad = None
if use_grad:
igrad = np.zeros(len(time)) * np.NaN
for t, loct in enumerate(temp):
slope, _, _, p_val, _ = stats.linregress(hgt.flatten(),
loct.flatten())
igrad[t] = slope if (p_val < 0.01) else np.NaN
gdir.write_monthly_climate_file(time, prcp[:, 1, 1], temp[:, 1, 1],
hgt[1, 1], ref_lon[1], ref_lat[1],
gradient=igrad)
# metadata
out = {'baseline_climate_source': source,
'baseline_hydro_yr_0': y0 + 1,
'baseline_hydro_yr_1': y1}
gdir.write_pickle(out, 'climate_info')
def mb_climate_on_height(gdir, heights, *, time_range=None, year_range=None):
"""Mass-balance climate of the glacier at a specific height
Reads the glacier's monthly climate data file and computes the
temperature "energies" (temp above 0) and solid precipitation at the
required height.
All MB parameters are considered here! (i.e. melt temp, precip scaling
factor, etc.)
Parameters
----------
gdir : GlacierDirectory
the glacier directory
heights: ndarray
a 1D array of the heights (in meter) where you want the data
time_range : [datetime, datetime], optional
default is to read all data but with this you
can provide a [t0, t1] bounds (inclusive).
year_range : [int, int], optional
Provide a [y0, y1] year range to get the data for specific
(hydrological) years only. Easier to use than the time bounds above.
Returns
-------
(time, tempformelt, prcpsol)::
- time: array of shape (nt,)
- tempformelt: array of shape (len(heights), nt)
- prcpsol: array of shape (len(heights), nt)
"""
if year_range is not None:
sm = cfg.PARAMS['hydro_month_' + gdir.hemisphere]
em = sm - 1 if (sm > 1) else 12
t0 = datetime.datetime(year_range[0]-1, sm, 1)
t1 = datetime.datetime(year_range[1], em, 1)
return mb_climate_on_height(gdir, heights, time_range=[t0, t1])
# Parameters
temp_all_solid = cfg.PARAMS['temp_all_solid']
temp_all_liq = cfg.PARAMS['temp_all_liq']
temp_melt = cfg.PARAMS['temp_melt']
prcp_fac = cfg.PARAMS['prcp_scaling_factor']
default_grad = cfg.PARAMS['temp_default_gradient']
g_minmax = cfg.PARAMS['temp_local_gradient_bounds']
# Read file
igrad = None
with utils.ncDataset(gdir.get_filepath('climate_monthly'), mode='r') as nc:
# time
time = nc.variables['time']
time = netCDF4.num2date(time[:], time.units)
if time_range is not None:
p0 = np.where(time == time_range[0])[0]
try:
p0 = p0[0]
except IndexError:
raise MassBalanceCalibrationError('time_range[0] not found in '
'file')
p1 = np.where(time == time_range[1])[0]
try:
p1 = p1[0]
except IndexError:
raise MassBalanceCalibrationError('time_range[1] not found in '
'file')
else:
p0 = 0
p1 = len(time)-1
time = time[p0:p1+1]
# Read timeseries
itemp = nc.variables['temp'][p0:p1+1]
iprcp = nc.variables['prcp'][p0:p1+1]
if 'gradient' in nc.variables:
igrad = nc.variables['gradient'][p0:p1+1]
# Security for stuff that can happen with local gradients
igrad = np.where(~np.isfinite(igrad), default_grad, igrad)
igrad = np.clip(igrad, g_minmax[0], g_minmax[1])
ref_hgt = nc.ref_hgt
# Default gradient?
if igrad is None:
igrad = itemp * 0 + default_grad
# Correct precipitation
iprcp *= prcp_fac
# For each height pixel:
# Compute temp and tempformelt (temperature above melting threshold)
npix = len(heights)
grad_temp = np.atleast_2d(igrad).repeat(npix, 0)
grad_temp *= (heights.repeat(len(time)).reshape(grad_temp.shape) - ref_hgt)
temp2d = np.atleast_2d(itemp).repeat(npix, 0) + grad_temp
temp2dformelt = temp2d - temp_melt
temp2dformelt = utils.clip_min(temp2dformelt, 0)
# Compute solid precipitation from total precipitation
prcpsol = np.atleast_2d(iprcp).repeat(npix, 0)
fac = 1 - (temp2d - temp_all_solid) / (temp_all_liq - temp_all_solid)
fac = np.clip(fac, 0, 1)
prcpsol = prcpsol * fac
return time, temp2dformelt, prcpsol
def mb_yearly_climate_on_height(gdir, heights, *,
year_range=None, flatten=False):
"""Yearly mass-balance climate of the glacier at a specific height
See also: mb_climate_on_height
Parameters
----------
gdir : GlacierDirectory
the glacier directory
heights: ndarray
a 1D array of the heights (in meter) where you want the data
year_range : [int, int], optional
Provide a [y0, y1] year range to get the data for specific
(hydrological) years only.
flatten : bool
for some applications (glacier average MB) it's ok to flatten the
data (average over height) prior to annual summing.
Returns
-------
(years, tempformelt, prcpsol)::
- years: array of shape (ny,)
- tempformelt: array of shape (len(heights), ny) (or ny if flatten
is set)
- prcpsol: array of shape (len(heights), ny) (or ny if flatten
is set)
"""
time, temp, prcp = mb_climate_on_height(gdir, heights,
year_range=year_range)
ny, r = divmod(len(time), 12)
if r != 0:
raise InvalidParamsError('Climate data should be N full years '
'exclusively')
# Last year gives the tone of the hydro year
years = np.arange(time[-1].year-ny+1, time[-1].year+1, 1)
if flatten:
# Spatial average
temp_yr = np.zeros(len(years))
prcp_yr = np.zeros(len(years))
temp = np.mean(temp, axis=0)
prcp = np.mean(prcp, axis=0)
for i, y in enumerate(years):
temp_yr[i] = np.sum(temp[i*12:(i+1)*12])
prcp_yr[i] = np.sum(prcp[i*12:(i+1)*12])
else:
# Annual prcp and temp for each point (no spatial average)
temp_yr = np.zeros((len(heights), len(years)))
prcp_yr = np.zeros((len(heights), len(years)))
for i, y in enumerate(years):
temp_yr[:, i] = np.sum(temp[:, i*12:(i+1)*12], axis=1)
prcp_yr[:, i] = np.sum(prcp[:, i*12:(i+1)*12], axis=1)
return years, temp_yr, prcp_yr
def mb_yearly_climate_on_glacier(gdir, *, year_range=None):
"""Yearly mass-balance climate at all glacier heights,
multiplied with the flowlines widths. (all in pix coords.)
See also: mb_climate_on_height
Parameters
----------
gdir : GlacierDirectory
the glacier directory
year_range : [int, int], optional
Provide a [y0, y1] year range to get the data for specific
(hydrological) years only.
Returns
-------
(years, tempformelt, prcpsol)::
- years: array of shape (ny)
- tempformelt: array of shape (ny)
- prcpsol: array of shape (ny)
"""
flowlines = gdir.read_pickle('inversion_flowlines')
heights = np.array([])
widths = np.array([])
for fl in flowlines:
heights = np.append(heights, fl.surface_h)
widths = np.append(widths, fl.widths)
years, temp, prcp = mb_yearly_climate_on_height(gdir, heights,
year_range=year_range,
flatten=False)
temp = np.average(temp, axis=0, weights=widths)
prcp = np.average(prcp, axis=0, weights=widths)
return years, temp, prcp
@entity_task(log, writes=['climate_info'])
def glacier_mu_candidates(gdir):
"""Computes the mu candidates, glacier wide.
For each 31 year-period centered on the year of interest, mu is is the
temperature sensitivity necessary for the glacier with its current shape
to be in equilibrium with its climate.
This task is just for documentation and testing! It is not used in
production anymore.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
the glacier directory to process
"""
warnings.warn('The task `glacier_mu_candidates` is deprecated. It should '
'only be used for testing.', DeprecationWarning)
mu_hp = int(cfg.PARAMS['mu_star_halfperiod'])
# Only get the years were we consider looking for tstar
y0, y1 = cfg.PARAMS['tstar_search_window']
ci = gdir.read_pickle('climate_info')
y0 = y0 or ci['baseline_hydro_yr_0']
y1 = y1 or ci['baseline_hydro_yr_1']
years, temp_yr, prcp_yr = mb_yearly_climate_on_glacier(gdir,
year_range=[y0, y1])
# Compute mu for each 31-yr climatological period
ny = len(years)
mu_yr_clim = np.zeros(ny) * np.NaN
for i, y in enumerate(years):
# Ignore begin and end
if ((i-mu_hp) < 0) or ((i+mu_hp) >= ny):
continue
t_avg = np.mean(temp_yr[i-mu_hp:i+mu_hp+1])
if t_avg > 1e-3: # if too cold no melt possible
prcp_ts = prcp_yr[i-mu_hp:i+mu_hp+1]
mu_yr_clim[i] = np.mean(prcp_ts) / t_avg
# Check that we found a least one mustar
if np.sum(np.isfinite(mu_yr_clim)) < 1:
raise MassBalanceCalibrationError('({}) no mustar candidates found.'
.format(gdir.rgi_id))
# Write
ci['mu_candidates_glacierwide'] = pd.Series(data=mu_yr_clim, index=years)
gdir.write_pickle(ci, 'climate_info')
@entity_task(log, writes=['climate_info'])
def t_star_from_refmb(gdir, mbdf=None, glacierwide=None,
write_diagnostics=False):
"""Computes the ref t* for the glacier, given a series of MB measurements.
Parameters
----------
gdir : oggm.GlacierDirectory
mbdf: a pd.Series containing the observed MB data indexed by year
if None, read automatically from the reference data
Returns
-------
A dict: {t_star:[], bias:[]}
"""
from oggm.core.massbalance import MultipleFlowlineMassBalance
if glacierwide is None:
glacierwide = cfg.PARAMS['tstar_search_glacierwide']
# Be sure we have no marine terminating glacier
assert not gdir.is_tidewater
# Reference time series
if mbdf is None:
mbdf = gdir.get_ref_mb_data()['ANNUAL_BALANCE']
# which years to look at
ref_years = mbdf.index.values
# Average oberved mass-balance
ref_mb = np.mean(mbdf)
# Compute one mu candidate per year and the associated statistics
# Only get the years were we consider looking for tstar
y0, y1 = cfg.PARAMS['tstar_search_window']
ci = gdir.read_pickle('climate_info')
y0 = y0 or ci['baseline_hydro_yr_0']
y1 = y1 or ci['baseline_hydro_yr_1']
years = np.arange(y0, y1+1)
ny = len(years)
mu_hp = int(cfg.PARAMS['mu_star_halfperiod'])
mb_per_mu = pd.Series(index=years)
if glacierwide:
# The old (but fast) method to find t*
_, temp, prcp = mb_yearly_climate_on_glacier(gdir, year_range=[y0, y1])
# which years to look at
selind = np.searchsorted(years, mbdf.index)
sel_temp = temp[selind]
sel_prcp = prcp[selind]
sel_temp = np.mean(sel_temp)
sel_prcp = np.mean(sel_prcp)
for i, y in enumerate(years):
# Ignore begin and end
if ((i - mu_hp) < 0) or ((i + mu_hp) >= ny):
continue
# Compute the mu candidate
t_avg = np.mean(temp[i - mu_hp:i + mu_hp + 1])
if t_avg < 1e-3: # if too cold no melt possible
continue
mu = np.mean(prcp[i - mu_hp:i + mu_hp + 1]) / t_avg
# Apply it
mb_per_mu[y] = np.mean(sel_prcp - mu * sel_temp)
else:
# The new (but slow) method to find t*
# Compute mu for each 31-yr climatological period
fls = gdir.read_pickle('inversion_flowlines')
for i, y in enumerate(years):
# Ignore begin and end
if ((i-mu_hp) < 0) or ((i+mu_hp) >= ny):
continue
# Calibrate the mu for this year
for fl in fls:
fl.mu_star_is_valid = False
try:
# TODO: this is slow and can be highly optimised
# it reads the same data over and over again
_recursive_mu_star_calibration(gdir, fls, y, first_call=True)
# Compute the MB with it
mb_mod = MultipleFlowlineMassBalance(gdir, fls, bias=0,
check_calib_params=False)
mb_ts = mb_mod.get_specific_mb(fls=fls, year=ref_years)
mb_per_mu[y] = np.mean(mb_ts)
except MassBalanceCalibrationError:
pass
# Diff to reference
diff = (mb_per_mu - ref_mb).dropna()
if len(diff) == 0:
raise MassBalanceCalibrationError('No single valid mu candidate for '
'this glacier!')
# Here we used to keep all possible mu* in order to later select
# them based on some distance search algorithms.
# (revision 81bc0923eab6301306184d26462f932b72b84117)
#
# As of Jul 2018, we will now stop this non-sense:
# out of all mu*, let's just pick the one with the smallest bias.
# It doesn't make much sense, but the same is true for other methods
# as well -> this is how Ben used to do it, and he is clever
# Another way would be to pick the closest to today or something
amin = np.abs(diff).idxmin()
# Write
d = gdir.read_pickle('climate_info')
d['t_star'] = amin
d['bias'] = diff[amin]
if write_diagnostics:
d['avg_mb_per_mu'] = mb_per_mu
d['avg_ref_mb'] = ref_mb
gdir.write_pickle(d, 'climate_info')
return {'t_star': amin, 'bias': diff[amin]}
def calving_mb(gdir):
"""Calving mass-loss in specific MB equivalent.
This is necessary to compute mu star.
"""
if not gdir.is_tidewater:
return 0.
# Ok. Just take the calving rate from cfg and change its units
# Original units: km3 a-1, to change to mm a-1 (units of specific MB)
rho = cfg.PARAMS['ice_density']
return gdir.inversion_calving_rate * 1e9 * rho / gdir.rgi_area_m2
def _fallback_local_t_star(gdir):
"""A Fallback function if climate.local_t_star raises an Error.
This function will still write a `local_mustar.json`, filled with NANs,
if climate.local_t_star fails and cfg.PARAMS['continue_on_error'] = True.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
the glacier directory to process
"""
# Scalars in a small dict for later
df = dict()
df['rgi_id'] = gdir.rgi_id