/
gis.py
1399 lines (1180 loc) · 49.2 KB
/
gis.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
""" Handling of the local glacier map and masks. Defines the first tasks
to be realized by any OGGM pre-processing workflow.
"""
# Built ins
import os
import logging
import warnings
from functools import partial
from distutils.version import LooseVersion
# External libs
import numpy as np
import shapely.ops
import pandas as pd
import shapely.geometry as shpg
import scipy.signal
from scipy.ndimage.measurements import label
from scipy.ndimage import binary_erosion
from scipy.ndimage.morphology import distance_transform_edt
from scipy.interpolate import griddata
from scipy import optimize as optimization
# Optional libs
try:
import salem
from salem.gis import transform_proj
except ImportError:
pass
try:
import pyproj
except ImportError:
pass
try:
import geopandas as gpd
except ImportError:
pass
try:
import skimage.draw as skdraw
except ImportError:
pass
try:
import rasterio
from rasterio.warp import reproject, Resampling
from rasterio.mask import mask as riomask
try:
# rasterio V > 1.0
from rasterio.merge import merge as merge_tool
except ImportError:
from rasterio.tools.merge import merge as merge_tool
except ImportError:
pass
# Locals
from oggm import entity_task, utils
import oggm.cfg as cfg
from oggm.exceptions import (InvalidParamsError, InvalidGeometryError,
InvalidDEMError, GeometryError)
from oggm.utils import (tuple2int, get_topo_file, is_dem_source_available,
nicenumber, ncDataset, tolist)
# Module logger
log = logging.getLogger(__name__)
# Needed later
label_struct = np.ones((3, 3))
def _parse_source_text():
fp = os.path.join(os.path.abspath(os.path.dirname(cfg.__file__)),
'data', 'dem_sources.txt')
out = dict()
cur_key = None
with open(fp, 'r') as fr:
this_text = []
for l in fr.readlines():
l = l.strip()
if l and (l[0] == '[' and l[-1] == ']'):
if cur_key:
out[cur_key] = '\n'.join(this_text)
this_text = []
cur_key = l.strip('[]')
continue
this_text.append(l)
out[cur_key] = '\n'.join(this_text)
return out
DEM_SOURCE_INFO = _parse_source_text()
def gaussian_blur(in_array, size):
"""Applies a Gaussian filter to a 2d array.
Parameters
----------
in_array : numpy.array
The array to smooth.
size : int
The half size of the smoothing window.
Returns
-------
a smoothed numpy.array
"""
# expand in_array to fit edge of kernel
padded_array = np.pad(in_array, size, 'symmetric')
# build kernel
x, y = np.mgrid[-size:size + 1, -size:size + 1]
g = np.exp(-(x**2 / float(size) + y**2 / float(size)))
g = (g / g.sum()).astype(in_array.dtype)
# do the Gaussian blur
return scipy.signal.fftconvolve(padded_array, g, mode='valid')
def multi_to_poly(geometry, gdir=None):
"""Sometimes an RGI geometry is a multipolygon: this should not happen.
Parameters
----------
geometry : shpg.Polygon or shpg.MultiPolygon
the geometry to check
gdir : GlacierDirectory, optional
for logging
Returns
-------
the corrected geometry
"""
# Log
rid = gdir.rgi_id + ': ' if gdir is not None else ''
if 'Multi' in geometry.type:
parts = np.array(geometry)
for p in parts:
assert p.type == 'Polygon'
areas = np.array([p.area for p in parts])
parts = parts[np.argsort(areas)][::-1]
areas = areas[np.argsort(areas)][::-1]
# First case (was RGIV4):
# let's assume that one poly is exterior and that
# the other polygons are in fact interiors
exterior = parts[0].exterior
interiors = []
was_interior = 0
for p in parts[1:]:
if parts[0].contains(p):
interiors.append(p.exterior)
was_interior += 1
if was_interior > 0:
# We are done here, good
geometry = shpg.Polygon(exterior, interiors)
else:
# This happens for bad geometries. We keep the largest
geometry = parts[0]
if np.any(areas[1:] > (areas[0] / 4)):
log.warning('Geometry {} lost quite a chunk.'.format(rid))
if geometry.type != 'Polygon':
raise InvalidGeometryError('Geometry {} is not a Polygon.'.format(rid))
return geometry
def _interp_polygon(polygon, dx):
"""Interpolates an irregular polygon to a regular step dx.
Interior geometries are also interpolated if they are longer then 3*dx,
otherwise they are ignored.
Parameters
----------
polygon: The shapely.geometry.Polygon instance to interpolate
dx : the step (float)
Returns
-------
an interpolated shapely.geometry.Polygon class instance.
"""
# remove last (duplex) point to build a LineString from the LinearRing
line = shpg.LineString(np.asarray(polygon.exterior.xy).T)
e_line = []
for distance in np.arange(0.0, line.length, dx):
e_line.append(*line.interpolate(distance).coords)
e_line = shpg.LinearRing(e_line)
i_lines = []
for ipoly in polygon.interiors:
line = shpg.LineString(np.asarray(ipoly.xy).T)
if line.length < 3*dx:
continue
i_points = []
for distance in np.arange(0.0, line.length, dx):
i_points.append(*line.interpolate(distance).coords)
i_lines.append(shpg.LinearRing(i_points))
return shpg.Polygon(e_line, i_lines)
def _polygon_to_pix(polygon):
"""Transforms polygon coordinates to integer pixel coordinates. It makes
the geometry easier to handle and reduces the number of points.
Parameters
----------
polygon: the shapely.geometry.Polygon instance to transform.
Returns
-------
a shapely.geometry.Polygon class instance.
"""
def project(x, y):
return np.rint(x).astype(np.int64), np.rint(y).astype(np.int64)
poly_pix = shapely.ops.transform(project, polygon)
# simple trick to correct invalid polys:
tmp = poly_pix.buffer(0)
# sometimes the glacier gets cut out in parts
if tmp.type == 'MultiPolygon':
# If only small arms are cut out, remove them
area = np.array([_tmp.area for _tmp in tmp])
_tokeep = np.argmax(area).item()
tmp = tmp[_tokeep]
# check that the other parts really are small,
# otherwise replace tmp with something better
area = area / area[_tokeep]
for _a in area:
if _a != 1 and _a > 0.05:
# these are extremely thin glaciers
# eg. RGI40-11.01381 RGI40-11.01697 params.d1 = 5. and d2 = 8.
# make them bigger until its ok
for b in np.arange(0., 1., 0.01):
tmp = shapely.ops.transform(project, polygon.buffer(b))
tmp = tmp.buffer(0)
if tmp.type == 'MultiPolygon':
continue
if tmp.is_valid:
break
if b == 0.99:
raise InvalidGeometryError('This glacier geometry is not '
'valid.')
if not tmp.is_valid:
raise InvalidGeometryError('This glacier geometry is not valid.')
return tmp
@entity_task(log, writes=['glacier_grid', 'dem', 'outlines'])
def define_glacier_region(gdir, entity=None):
"""Very first task: define the glacier's local grid.
Defines the local projection (Transverse Mercator), centered on the
glacier. There is some options to set the resolution of the local grid.
It can be adapted depending on the size of the glacier with::
dx (m) = d1 * AREA (km) + d2 ; clipped to dmax
or be set to a fixed value. See ``params.cfg`` for setting these options.
Default values of the adapted mode lead to a resolution of 50 m for
Hintereisferner, which is approx. 8 km2 large.
After defining the grid, the topography and the outlines of the glacier
are transformed into the local projection. The default interpolation for
the topography is `cubic`.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
where to write the data
entity : geopandas.GeoSeries
the glacier geometry to process
"""
# Make a local glacier map
proj_params = dict(name='tmerc', lat_0=0., lon_0=gdir.cenlon,
k=0.9996, x_0=0, y_0=0, datum='WGS84')
proj4_str = "+proj={name} +lat_0={lat_0} +lon_0={lon_0} +k={k} " \
"+x_0={x_0} +y_0={y_0} +datum={datum}".format(**proj_params)
proj_in = pyproj.Proj("+init=EPSG:4326", preserve_units=True)
proj_out = pyproj.Proj(proj4_str, preserve_units=True)
project = partial(transform_proj, proj_in, proj_out)
# transform geometry to map
geometry = shapely.ops.transform(project, entity['geometry'])
geometry = multi_to_poly(geometry, gdir=gdir)
xx, yy = geometry.exterior.xy
# Save transformed geometry to disk
entity = entity.copy()
entity['geometry'] = geometry
# Avoid fiona bug: https://github.com/Toblerity/Fiona/issues/365
for k, s in entity.iteritems():
if type(s) in [np.int32, np.int64]:
entity[k] = int(s)
towrite = gpd.GeoDataFrame(entity).T
towrite.crs = proj4_str
# Delete the source before writing
if 'DEM_SOURCE' in towrite:
del towrite['DEM_SOURCE']
# Define glacier area to use
area = entity['Area']
# Do we want to use the RGI area or ours?
if not cfg.PARAMS['use_rgi_area']:
# Update Area
area = geometry.area * 1e-6
entity['Area'] = area
towrite['Area'] = area
# Write shapefile
gdir.write_shapefile(towrite, 'outlines')
# Also transform the intersects if necessary
gdf = cfg.PARAMS['intersects_gdf']
if len(gdf) > 0:
gdf = gdf.loc[((gdf.RGIId_1 == gdir.rgi_id) |
(gdf.RGIId_2 == gdir.rgi_id))]
if len(gdf) > 0:
gdf = salem.transform_geopandas(gdf, to_crs=proj_out)
if hasattr(gdf.crs, 'srs'):
# salem uses pyproj
gdf.crs = gdf.crs.srs
gdir.write_shapefile(gdf, 'intersects')
else:
# Sanity check
if cfg.PARAMS['use_intersects']:
raise InvalidParamsError('You seem to have forgotten to set the '
'intersects file for this run. OGGM '
'works better with such a file. If you '
'know what your are doing, set '
"cfg.PARAMS['use_intersects'] = False to "
"suppress this error.")
# 6. choose a spatial resolution with respect to the glacier area
dxmethod = cfg.PARAMS['grid_dx_method']
if dxmethod == 'linear':
dx = np.rint(cfg.PARAMS['d1'] * area + cfg.PARAMS['d2'])
elif dxmethod == 'square':
dx = np.rint(cfg.PARAMS['d1'] * np.sqrt(area) + cfg.PARAMS['d2'])
elif dxmethod == 'fixed':
dx = np.rint(cfg.PARAMS['fixed_dx'])
else:
raise InvalidParamsError('grid_dx_method not supported: {}'
.format(dxmethod))
# Additional trick for varying dx
if dxmethod in ['linear', 'square']:
dx = utils.clip_scalar(dx, cfg.PARAMS['d2'], cfg.PARAMS['dmax'])
log.debug('(%s) area %.2f km, dx=%.1f', gdir.rgi_id, area, dx)
# Safety check
border = cfg.PARAMS['border']
if border > 1000:
raise InvalidParamsError("You have set a cfg.PARAMS['border'] value "
"of {}. ".format(cfg.PARAMS['border']) +
'This a very large value, which is '
'currently not supported in OGGM.')
# For tidewater glaciers we force border to 10
if gdir.is_tidewater and cfg.PARAMS['clip_tidewater_border']:
border = 10
# Corners, incl. a buffer of N pix
ulx = np.min(xx) - border * dx
lrx = np.max(xx) + border * dx
uly = np.max(yy) + border * dx
lry = np.min(yy) - border * dx
# n pixels
nx = np.int((lrx - ulx) / dx)
ny = np.int((uly - lry) / dx)
# Back to lon, lat for DEM download/preparation
tmp_grid = salem.Grid(proj=proj_out, nxny=(nx, ny), x0y0=(ulx, uly),
dxdy=(dx, -dx), pixel_ref='corner')
minlon, maxlon, minlat, maxlat = tmp_grid.extent_in_crs(crs=salem.wgs84)
# Open DEM
source = entity.DEM_SOURCE if hasattr(entity, 'DEM_SOURCE') else None
if not is_dem_source_available(source,
(minlon, maxlon),
(minlat, maxlat)):
raise InvalidDEMError('Source: {} not available for glacier {}'
.format(source, gdir.rgi_id))
dem_list, dem_source = get_topo_file((minlon, maxlon), (minlat, maxlat),
rgi_region=gdir.rgi_region,
rgi_subregion=gdir.rgi_subregion,
dx_meter=dx,
source=source)
log.debug('(%s) DEM source: %s', gdir.rgi_id, dem_source)
log.debug('(%s) N DEM Files: %s', gdir.rgi_id, len(dem_list))
# A glacier area can cover more than one tile:
if len(dem_list) == 1:
dem_dss = [rasterio.open(dem_list[0])] # if one tile, just open it
dem_data = rasterio.band(dem_dss[0], 1)
if LooseVersion(rasterio.__version__) >= LooseVersion('1.0'):
src_transform = dem_dss[0].transform
else:
src_transform = dem_dss[0].affine
else:
dem_dss = [rasterio.open(s) for s in dem_list] # list of rasters
dem_data, src_transform = merge_tool(dem_dss) # merged rasters
# Use Grid properties to create a transform (see rasterio cookbook)
dst_transform = rasterio.transform.from_origin(
ulx, uly, dx, dx # sign change (2nd dx) is done by rasterio.transform
)
# Set up profile for writing output
profile = dem_dss[0].profile
profile.update({
'crs': proj4_str,
'transform': dst_transform,
'width': nx,
'height': ny
})
# Could be extended so that the cfg file takes all Resampling.* methods
if cfg.PARAMS['topo_interp'] == 'bilinear':
resampling = Resampling.bilinear
elif cfg.PARAMS['topo_interp'] == 'cubic':
resampling = Resampling.cubic
else:
raise InvalidParamsError('{} interpolation not understood'
.format(cfg.PARAMS['topo_interp']))
dem_reproj = gdir.get_filepath('dem')
profile.pop('blockxsize', None)
profile.pop('blockysize', None)
profile.pop('compress', None)
nodata = dem_dss[0].meta.get('nodata', None)
if source == 'TANDEM' and nodata is None:
# badly tagged geotiffs, let's do it ourselves
nodata = -32767
with rasterio.open(dem_reproj, 'w', **profile) as dest:
dst_array = np.empty((ny, nx), dtype=dem_dss[0].dtypes[0])
reproject(
# Source parameters
source=dem_data,
src_crs=dem_dss[0].crs,
src_transform=src_transform,
src_nodata=nodata,
# Destination parameters
destination=dst_array,
dst_transform=dst_transform,
dst_crs=proj4_str,
dst_nodata=nodata,
# Configuration
resampling=resampling)
dest.write(dst_array, 1)
for dem_ds in dem_dss:
dem_ds.close()
# Glacier grid
x0y0 = (ulx+dx/2, uly-dx/2) # To pixel center coordinates
glacier_grid = salem.Grid(proj=proj_out, nxny=(nx, ny), dxdy=(dx, -dx),
x0y0=x0y0)
glacier_grid.to_json(gdir.get_filepath('glacier_grid'))
# Write DEM source info
gdir.add_to_diagnostics('dem_source', dem_source)
source_txt = DEM_SOURCE_INFO.get(dem_source, dem_source)
with open(gdir.get_filepath('dem_source'), 'w') as fw:
fw.write(source_txt)
fw.write('\n\n')
fw.write('# Data files\n\n')
for fname in dem_list:
fw.write('{}\n'.format(os.path.basename(fname)))
@entity_task(log, writes=['gridded_data', 'geometries'])
def glacier_masks(gdir):
"""Makes a gridded mask of the glacier outlines and topography.
This function fills holes in the source DEM and produces smoothed gridded
topography and glacier outline arrays. These are the ones which will later
be used to determine bed and surface height.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
where to write the data
"""
# In case nominal, just raise
if gdir.is_nominal:
raise GeometryError('{} is a nominal glacier.'.format(gdir.rgi_id))
# open srtm tif-file:
dem_dr = rasterio.open(gdir.get_filepath('dem'), 'r', driver='GTiff')
dem = dem_dr.read(1).astype(rasterio.float32)
# Grid
nx = dem_dr.width
ny = dem_dr.height
assert nx == gdir.grid.nx
assert ny == gdir.grid.ny
# Correct the DEM
# Currently we just do a linear interp -- filling is totally shit anyway
min_z = -999.
dem[dem <= min_z] = np.NaN
isfinite = np.isfinite(dem)
if np.all(~isfinite):
raise InvalidDEMError('Not a single valid grid point in DEM')
if np.any(~isfinite):
xx, yy = gdir.grid.ij_coordinates
pnan = np.nonzero(~isfinite)
pok = np.nonzero(isfinite)
points = np.array((np.ravel(yy[pok]), np.ravel(xx[pok]))).T
inter = np.array((np.ravel(yy[pnan]), np.ravel(xx[pnan]))).T
try:
dem[pnan] = griddata(points, np.ravel(dem[pok]), inter,
method='linear')
except ValueError:
raise InvalidDEMError('DEM interpolation not possible.')
log.warning(gdir.rgi_id + ': DEM needed interpolation.')
gdir.add_to_diagnostics('dem_needed_interpolation', True)
gdir.add_to_diagnostics('dem_invalid_perc', len(pnan[0]) / (nx*ny))
isfinite = np.isfinite(dem)
if np.any(~isfinite):
# this happens when extrapolation is needed
# see how many percent of the dem
if np.sum(~isfinite) > (0.5 * nx * ny):
log.warning('({}) many NaNs in DEM'.format(gdir.rgi_id))
xx, yy = gdir.grid.ij_coordinates
pnan = np.nonzero(~isfinite)
pok = np.nonzero(isfinite)
points = np.array((np.ravel(yy[pok]), np.ravel(xx[pok]))).T
inter = np.array((np.ravel(yy[pnan]), np.ravel(xx[pnan]))).T
try:
dem[pnan] = griddata(points, np.ravel(dem[pok]), inter,
method='nearest')
except ValueError:
raise InvalidDEMError('DEM extrapolation not possible.')
log.warning(gdir.rgi_id + ': DEM needed extrapolation.')
gdir.add_to_diagnostics('dem_needed_extrapolation', True)
gdir.add_to_diagnostics('dem_extrapol_perc', len(pnan[0]) / (nx*ny))
if np.min(dem) == np.max(dem):
raise InvalidDEMError('({}) min equal max in the DEM.'
.format(gdir.rgi_id))
# Projection
if LooseVersion(rasterio.__version__) >= LooseVersion('1.0'):
transf = dem_dr.transform
else:
transf = dem_dr.affine
x0 = transf[2] # UL corner
y0 = transf[5] # UL corner
dx = transf[0]
dy = transf[4] # Negative
if not (np.allclose(dx, -dy) or np.allclose(dx, gdir.grid.dx) or
np.allclose(y0, gdir.grid.corner_grid.y0, atol=1e-2) or
np.allclose(x0, gdir.grid.corner_grid.x0, atol=1e-2)):
raise InvalidDEMError('DEM file and Salem Grid do not match!')
dem_dr.close()
# Clip topography to 0 m a.s.l.
utils.clip_min(dem, 0, out=dem)
# Smooth DEM?
if cfg.PARAMS['smooth_window'] > 0.:
gsize = np.rint(cfg.PARAMS['smooth_window'] / dx)
smoothed_dem = gaussian_blur(dem, np.int(gsize))
else:
smoothed_dem = dem.copy()
if not np.all(np.isfinite(smoothed_dem)):
raise InvalidDEMError('({}) NaN in smoothed DEM'.format(gdir.rgi_id))
# Geometries
geometry = gdir.read_shapefile('outlines').geometry[0]
# Interpolate shape to a regular path
glacier_poly_hr = _interp_polygon(geometry, gdir.grid.dx)
# Transform geometry into grid coordinates
# It has to be in pix center coordinates because of how skimage works
def proj(x, y):
grid = gdir.grid.center_grid
return grid.transform(x, y, crs=grid.proj)
glacier_poly_hr = shapely.ops.transform(proj, glacier_poly_hr)
# simple trick to correct invalid polys:
# http://stackoverflow.com/questions/20833344/
# fix-invalid-polygon-python-shapely
glacier_poly_hr = glacier_poly_hr.buffer(0)
if not glacier_poly_hr.is_valid:
raise InvalidGeometryError('This glacier geometry is not valid.')
# Rounded nearest pix
glacier_poly_pix = _polygon_to_pix(glacier_poly_hr)
if glacier_poly_pix.exterior is None:
raise InvalidGeometryError('Problem in converting glacier geometry '
'to grid resolution.')
# Compute the glacier mask (currently: center pixels + touched)
nx, ny = gdir.grid.nx, gdir.grid.ny
glacier_mask = np.zeros((ny, nx), dtype=np.uint8)
glacier_ext = np.zeros((ny, nx), dtype=np.uint8)
(x, y) = glacier_poly_pix.exterior.xy
glacier_mask[skdraw.polygon(np.array(y), np.array(x))] = 1
for gint in glacier_poly_pix.interiors:
x, y = tuple2int(gint.xy)
glacier_mask[skdraw.polygon(y, x)] = 0
glacier_mask[y, x] = 0 # on the nunataks, no
x, y = tuple2int(glacier_poly_pix.exterior.xy)
glacier_mask[y, x] = 1
glacier_ext[y, x] = 1
# Because of the 0 values at nunataks boundaries, some "Ice Islands"
# can happen within nunataks (e.g.: RGI40-11.00062)
# See if we can filter them out easily
regions, nregions = label(glacier_mask, structure=label_struct)
if nregions > 1:
log.debug('(%s) we had to cut an island in the mask', gdir.rgi_id)
# Check the size of those
region_sizes = [np.sum(regions == r) for r in np.arange(1, nregions+1)]
am = np.argmax(region_sizes)
# Check not a strange glacier
sr = region_sizes.pop(am)
for ss in region_sizes:
assert (ss / sr) < 0.1
glacier_mask[:] = 0
glacier_mask[np.where(regions == (am+1))] = 1
# Last sanity check based on the masked dem
tmp_max = np.max(dem[np.where(glacier_mask == 1)])
tmp_min = np.min(dem[np.where(glacier_mask == 1)])
if tmp_max < (tmp_min + 1):
raise InvalidDEMError('({}) min equal max in the masked DEM.'
.format(gdir.rgi_id))
# write out the grids in the netcdf file
nc = gdir.create_gridded_ncdf_file('gridded_data')
v = nc.createVariable('topo', 'f4', ('y', 'x', ), zlib=True)
v.units = 'm'
v.long_name = 'DEM topography'
v[:] = dem
v = nc.createVariable('topo_smoothed', 'f4', ('y', 'x', ), zlib=True)
v.units = 'm'
v.long_name = ('DEM topography smoothed '
'with radius: {:.1} m'.format(cfg.PARAMS['smooth_window']))
v[:] = smoothed_dem
v = nc.createVariable('glacier_mask', 'i1', ('y', 'x', ), zlib=True)
v.units = '-'
v.long_name = 'Glacier mask'
v[:] = glacier_mask
v = nc.createVariable('glacier_ext', 'i1', ('y', 'x', ), zlib=True)
v.units = '-'
v.long_name = 'Glacier external boundaries'
v[:] = glacier_ext
# add some meta stats and close
nc.max_h_dem = np.max(dem)
nc.min_h_dem = np.min(dem)
dem_on_g = dem[np.where(glacier_mask)]
nc.max_h_glacier = np.max(dem_on_g)
nc.min_h_glacier = np.min(dem_on_g)
nc.close()
geometries = dict()
geometries['polygon_hr'] = glacier_poly_hr
geometries['polygon_pix'] = glacier_poly_pix
geometries['polygon_area'] = geometry.area
gdir.write_pickle(geometries, 'geometries')
@entity_task(log, writes=['gridded_data', 'hypsometry'])
def simple_glacier_masks(gdir):
"""Compute glacier masks based on much simpler rules than OGGM's default.
This is therefore more robust: we use this function to compute glacier
hypsometries.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
where to write the data
"""
# open srtm tif-file:
dem_dr = rasterio.open(gdir.get_filepath('dem'), 'r', driver='GTiff')
dem = dem_dr.read(1).astype(rasterio.float32)
# Grid
nx = dem_dr.width
ny = dem_dr.height
assert nx == gdir.grid.nx
assert ny == gdir.grid.ny
# Correct the DEM
# Currently we just do a linear interp -- filling is totally shit anyway
min_z = -999.
dem[dem <= min_z] = np.NaN
isfinite = np.isfinite(dem)
if np.all(~isfinite):
raise InvalidDEMError('Not a single valid grid point in DEM')
if np.any(~isfinite):
xx, yy = gdir.grid.ij_coordinates
pnan = np.nonzero(~isfinite)
pok = np.nonzero(isfinite)
points = np.array((np.ravel(yy[pok]), np.ravel(xx[pok]))).T
inter = np.array((np.ravel(yy[pnan]), np.ravel(xx[pnan]))).T
try:
dem[pnan] = griddata(points, np.ravel(dem[pok]), inter,
method='linear')
except ValueError:
raise InvalidDEMError('DEM interpolation not possible.')
log.warning(gdir.rgi_id + ': DEM needed interpolation.')
gdir.add_to_diagnostics('dem_needed_interpolation', True)
gdir.add_to_diagnostics('dem_invalid_perc', len(pnan[0]) / (nx*ny))
isfinite = np.isfinite(dem)
if np.any(~isfinite):
# this happens when extrapolation is needed
# see how many percent of the dem
if np.sum(~isfinite) > (0.5 * nx * ny):
log.warning('({}) many NaNs in DEM'.format(gdir.rgi_id))
xx, yy = gdir.grid.ij_coordinates
pnan = np.nonzero(~isfinite)
pok = np.nonzero(isfinite)
points = np.array((np.ravel(yy[pok]), np.ravel(xx[pok]))).T
inter = np.array((np.ravel(yy[pnan]), np.ravel(xx[pnan]))).T
try:
dem[pnan] = griddata(points, np.ravel(dem[pok]), inter,
method='nearest')
except ValueError:
raise InvalidDEMError('DEM extrapolation not possible.')
log.warning(gdir.rgi_id + ': DEM needed extrapolation.')
gdir.add_to_diagnostics('dem_needed_extrapolation', True)
gdir.add_to_diagnostics('dem_extrapol_perc', len(pnan[0]) / (nx*ny))
if np.min(dem) == np.max(dem):
raise InvalidDEMError('({}) min equal max in the DEM.'
.format(gdir.rgi_id))
# Proj
if LooseVersion(rasterio.__version__) >= LooseVersion('1.0'):
transf = dem_dr.transform
else:
raise ImportError('This task needs rasterio >= 1.0 to work properly')
x0 = transf[2] # UL corner
y0 = transf[5] # UL corner
dx = transf[0]
dy = transf[4] # Negative
assert dx == -dy
assert dx == gdir.grid.dx
assert y0 == gdir.grid.corner_grid.y0
assert x0 == gdir.grid.corner_grid.x0
profile = dem_dr.profile
dem_dr.close()
# Clip topography to 0 m a.s.l.
utils.clip_min(dem, 0, out=dem)
# Smooth DEM?
if cfg.PARAMS['smooth_window'] > 0.:
gsize = np.rint(cfg.PARAMS['smooth_window'] / dx)
smoothed_dem = gaussian_blur(dem, np.int(gsize))
else:
smoothed_dem = dem.copy()
if not np.all(np.isfinite(smoothed_dem)):
raise InvalidDEMError('({}) NaN in smoothed DEM'.format(gdir.rgi_id))
# Geometries
geometry = gdir.read_shapefile('outlines').geometry[0]
# simple trick to correct invalid polys:
# http://stackoverflow.com/questions/20833344/
# fix-invalid-polygon-python-shapely
geometry = geometry.buffer(0)
if not geometry.is_valid:
raise InvalidDEMError('This glacier geometry is not valid.')
# Compute the glacier mask using rasterio
# Small detour as mask only accepts DataReader objects
with rasterio.io.MemoryFile() as memfile:
with memfile.open(**profile) as dataset:
dataset.write(dem.astype(np.int16)[np.newaxis, ...])
dem_data = rasterio.open(memfile.name)
masked_dem, _ = riomask(dem_data, [shpg.mapping(geometry)],
filled=False)
glacier_mask = ~masked_dem[0, ...].mask
# Smame without nunataks
with rasterio.io.MemoryFile() as memfile:
with memfile.open(**profile) as dataset:
dataset.write(dem.astype(np.int16)[np.newaxis, ...])
dem_data = rasterio.open(memfile.name)
poly = shpg.mapping(shpg.Polygon(geometry.exterior))
masked_dem, _ = riomask(dem_data, [poly],
filled=False)
glacier_mask_nonuna = ~masked_dem[0, ...].mask
# Glacier exterior excluding nunataks
erode = binary_erosion(glacier_mask_nonuna)
glacier_ext = glacier_mask_nonuna ^ erode
glacier_ext = np.where(glacier_mask_nonuna, glacier_ext, 0)
# Last sanity check based on the masked dem
tmp_max = np.max(dem[glacier_mask])
tmp_min = np.min(dem[glacier_mask])
if tmp_max < (tmp_min + 1):
raise InvalidDEMError('({}) min equal max in the masked DEM.'
.format(gdir.rgi_id))
# hypsometry
bsize = 50.
dem_on_ice = dem[glacier_mask]
bins = np.arange(nicenumber(dem_on_ice.min(), bsize, lower=True),
nicenumber(dem_on_ice.max(), bsize) + 0.01, bsize)
h, _ = np.histogram(dem_on_ice, bins)
h = h / np.sum(h) * 1000 # in permil
# We want to convert the bins to ints but preserve their sum to 1000
# Start with everything rounded down, then round up the numbers with the
# highest fractional parts until the desired sum is reached.
hi = np.floor(h).astype(np.int)
hup = np.ceil(h).astype(np.int)
aso = np.argsort(hup - h)
for i in aso:
hi[i] = hup[i]
if np.sum(hi) == 1000:
break
# slope
sy, sx = np.gradient(dem, dx)
aspect = np.arctan2(np.mean(-sx[glacier_mask]), np.mean(sy[glacier_mask]))
aspect = np.rad2deg(aspect)
if aspect < 0:
aspect += 360
slope = np.arctan(np.sqrt(sx ** 2 + sy ** 2))
avg_slope = np.rad2deg(np.mean(slope[glacier_mask]))
# write
df = pd.DataFrame()
df['RGIId'] = [gdir.rgi_id]
df['GLIMSId'] = [gdir.glims_id]
df['Zmin'] = [dem_on_ice.min()]
df['Zmax'] = [dem_on_ice.max()]
df['Zmed'] = [np.median(dem_on_ice)]
df['Area'] = [gdir.rgi_area_km2]
df['Slope'] = [avg_slope]
df['Aspect'] = [aspect]
for b, bs in zip(hi, (bins[1:] + bins[:-1])/2):
df['{}'.format(np.round(bs).astype(int))] = [b]
df.to_csv(gdir.get_filepath('hypsometry'), index=False)
# write out the grids in the netcdf file
nc = gdir.create_gridded_ncdf_file('gridded_data')
v = nc.createVariable('topo', 'f4', ('y', 'x', ), zlib=True)
v.units = 'm'
v.long_name = 'DEM topography'
v[:] = dem
v = nc.createVariable('topo_smoothed', 'f4', ('y', 'x', ), zlib=True)
v.units = 'm'
v.long_name = ('DEM topography smoothed '
'with radius: {:.1} m'.format(cfg.PARAMS['smooth_window']))
v[:] = smoothed_dem
v = nc.createVariable('glacier_mask', 'i1', ('y', 'x', ), zlib=True)
v.units = '-'
v.long_name = 'Glacier mask'
v[:] = glacier_mask
v = nc.createVariable('glacier_ext', 'i1', ('y', 'x', ), zlib=True)
v.units = '-'
v.long_name = 'Glacier external boundaries'
v[:] = glacier_ext
# add some meta stats and close
nc.max_h_dem = np.max(dem)
nc.min_h_dem = np.min(dem)
dem_on_g = dem[np.where(glacier_mask)]
nc.max_h_glacier = np.max(dem_on_g)
nc.min_h_glacier = np.min(dem_on_g)
nc.close()
@entity_task(log, writes=['gridded_data'])
def gridded_attributes(gdir):
"""Adds attributes to the gridded file, useful for thickness interpolation.
This could be useful for distributed ice thickness models.
The raster data are added to the gridded_data file.
Parameters
----------
gdir : :py:class:`oggm.GlacierDirectory`
where to write the data
"""
# Variables
grids_file = gdir.get_filepath('gridded_data')
with ncDataset(grids_file) as nc:
topo_smoothed = nc.variables['topo_smoothed'][:]
glacier_mask = nc.variables['glacier_mask'][:]
# Glacier exterior including nunataks
erode = binary_erosion(glacier_mask)
glacier_ext = glacier_mask ^ erode
glacier_ext = np.where(glacier_mask == 1, glacier_ext, 0)
# Intersects between glaciers
gdfi = gpd.GeoDataFrame(columns=['geometry'])
if gdir.has_file('intersects'):
# read and transform to grid
gdf = gdir.read_shapefile('intersects')
salem.transform_geopandas(gdf, gdir.grid, inplace=True)
gdfi = pd.concat([gdfi, gdf[['geometry']]])
# Ice divide mask
# Probably not the fastest way to do this, but it works
dist = np.array([])
jj, ii = np.where(glacier_ext)
for j, i in zip(jj, ii):
dist = np.append(dist, np.min(gdfi.distance(shpg.Point(i, j))))
with warnings.catch_warnings():
warnings.filterwarnings("ignore", category=RuntimeWarning)
pok = np.where(dist <= 1)
glacier_ext_intersect = glacier_ext * 0
glacier_ext_intersect[jj[pok], ii[pok]] = 1
# Distance from border mask - Scipy does the job
dx = gdir.grid.dx
dis_from_border = 1 + glacier_ext_intersect - glacier_ext
dis_from_border = distance_transform_edt(dis_from_border) * dx
# Slope
glen_n = cfg.PARAMS['glen_n']
sy, sx = np.gradient(topo_smoothed, dx, dx)
slope = np.arctan(np.sqrt(sy**2 + sx**2))
slope_factor = utils.clip_scalar(slope,
np.deg2rad(cfg.PARAMS['min_slope']*4),
np.pi/2)
slope_factor = 1 / slope_factor**(glen_n / (glen_n+2))
aspect = np.arctan2(-sx, sy)
aspect[aspect < 0] += 2 * np.pi
with ncDataset(grids_file, 'a') as nc:
vn = 'glacier_ext_erosion'
if vn in nc.variables:
v = nc.variables[vn]
else:
v = nc.createVariable(vn, 'i1', ('y', 'x', ))
v.units = '-'
v.long_name = 'Glacier exterior with binary erosion method'
v[:] = glacier_ext
vn = 'ice_divides'
if vn in nc.variables:
v = nc.variables[vn]
else:
v = nc.createVariable(vn, 'i1', ('y', 'x', ))
v.units = '-'
v.long_name = 'Glacier ice divides'
v[:] = glacier_ext_intersect
vn = 'slope'
if vn in nc.variables:
v = nc.variables[vn]
else:
v = nc.createVariable(vn, 'f4', ('y', 'x', ))
v.units = 'rad'
v.long_name = 'Local slope based on smoothed topography'
v[:] = slope
vn = 'aspect'
if vn in nc.variables:
v = nc.variables[vn]