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
def multipolygon_to_polygon(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.info('Geometry {} lost quite a chunk.'.format(rid)) if geometry.type != 'Polygon': raise InvalidGeometryError('Geometry {} is not a Polygon.'.format(rid)) return geometry
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 """ # 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. dem = dem.clip(0) # 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')
def glacier_masks(gdir): """Makes a gridded mask of the glacier outlines that can be used by OGGM. For a more robust solution (not OGGM compatible) see simple_glacier_masks. 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)) if not os.path.exists(gdir.get_filepath('gridded_data')): # In a possible future, we might actually want to raise a # deprecation warning here process_dem(gdir) # 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 # Write geometries geometries = dict() geometries['polygon_hr'] = glacier_poly_hr geometries['polygon_pix'] = glacier_poly_pix geometries['polygon_area'] = geometry.area gdir.write_pickle(geometries, 'geometries') # write out the grids in the netcdf file with GriddedNcdfFile(gdir) as nc: if 'glacier_mask' not in nc.variables: v = nc.createVariable('glacier_mask', 'i1', ('y', 'x', ), zlib=True) v.units = '-' v.long_name = 'Glacier mask' else: v = nc.variables['glacier_mask'] v[:] = glacier_mask if 'glacier_ext' not in nc.variables: v = nc.createVariable('glacier_ext', 'i1', ('y', 'x', ), zlib=True) v.units = '-' v.long_name = 'Glacier external boundaries' else: v = nc.variables['glacier_ext'] v[:] = glacier_ext dem = nc.variables['topo'][:] valid_mask = nc.variables['topo_valid_mask'][:] # 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)) # Log DEM that needed processing within the glacier mask if gdir.get_diagnostics().get('dem_needed_interpolation', False): pnan = (valid_mask == 0) & glacier_mask gdir.add_to_diagnostics('dem_invalid_perc_in_mask', np.sum(pnan) / np.sum(glacier_mask)) # add some meta stats and close 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)