Prepare flow directions and related data from a DEM#
With HydroMT-Wflow, a user can choose to build a model in a geographic or projected coordinate system from an input Digital Elevation Model (DEM) and Flow Direction (flwdir) dataset.
While DEM data are often available, this is not the always the case for the flow directions (flwdir). We made the choice to build a Wflow model directly from user provided DEM and flwdir datasets rather than reprojecting a DEM and/or deriving flwdir on the fly. This is because deriving flow directions is often an iterative process to be sure the flow directions matche the terrain and river network. Note that for the best results the flwdir data should be derived from a high-res DEM (<100 m spatial resolution). The HydroMT-Wflow model builder will automatically resample the flwdir data to the model resolution.
Because of this, we prefer to provide this notebook as a possible pre-processing step before calling a build a Wflow model with HydroMt-Wflow. Here we use the different flow directions methods from HydroMT and PyFlwDir
Load dependencies#
[1]:
import os
import xarray as xr
import numpy as np
import geopandas as gpd
# pyflwdir
import pyflwdir
# hydromt
from hydromt import DataCatalog, flw
# plot
import matplotlib.pyplot as plt
from matplotlib import cm, colors
Deriving flow directions from Elevation data#
In this example we will use the merit_hydro
data in the pre-defined artifact_data
catalog of HydroMT.
The typical starting point to derive flow direction is a DEM and a river network vector file. However, just for this example were we already have flow direction data, we will derive the river network from the flow direction data. And then use the DEM and the river network to derive flow directions.
[2]:
# read data
data_catalog = DataCatalog("artifact_data")
ds_hydro_org = data_catalog.get_rasterdataset(
"merit_hydro", variables=["elevtn", "flwdir", "uparea", "strord"], bbox=[12, 46.0, 12.3, 46.2]
)
# derive river network
# Note: this is typically not needed, instead a user supplied river network is used
flwdir_org = flw.flwdir_from_da(ds_hydro_org["flwdir"], ftype="infer", check_ftype=True)
gdf_riv_org = gpd.GeoDataFrame.from_features(
flwdir_org.streams(uparea=ds_hydro_org["uparea"].values, min_sto=3),
crs=ds_hydro_org.raster.crs
)
To derive flow directions from a DEM, you can use the hydromt.flw.d8_from_dem method of HydroMT.
This method derives D8 flow directions grid from an elevation grid and allows several options to the users: - river burning: provide a river vector layer gdf_riv
with uparea
(km2) or rivdph
column to guide the derivation of flow directions. This is important to get the correct flow directions, especially in flat areas. - outlets: outlets can be defined at edge
s of the grid (defualt) or force all flow to go to the minimum elevation point min
. The latter only makes
sense if your DEM only is masked to the catchment. Additionnally, the user can also force specific pits locations via idxs_pit
. - depression filling: local depressions are filled based on their lowest pour point level if the pour point depth is smaller than the maximum pour point depth max_depth
, otherwise the lowest elevation in the depression becomes a pit. By default max_depth
is set to -1 m filling all local depressions.
[3]:
# Derive flow directions with outlets at the edges -> this is the default
da_flwdir = flw.d8_from_dem(
da_elv=ds_hydro_org["elevtn"],
max_depth=-1, # max depression poir point depth; -1 means no local pits
outlets="edge", # option: "edge" (default), "min", "idxs_pit"
idxs_pit=None,
gdf_riv=gdf_riv_org, # user supplied river network to aid flow direction derivation
riv_burn_method="uparea", # options: "fixed" (default), "rivdph", "uparea"
# riv_depth=5, # fixed river depth in meters, only used if riv_burn_method="fixed"
# **kwargs to be passed to pyflwdir.dem.fill_depressions
)
# convert to vector for plotting based no minimum stream order
flwdir = flw.flwdir_from_da(
da_flwdir, ftype="infer", check_ftype=True
)
gdf_riv = gpd.GeoDataFrame.from_features(
flwdir.streams(min_sto=3), crs=da_flwdir.raster.crs
)
Let’s compare the original with the derived flow direction data. Note that we do not see a big difference here, because the river network file that we used was created from the original flow directions.
[4]:
# plot
# initialize image with geoaxes
fig = plt.figure(figsize=(8, 8))
ax = fig.add_subplot()
# set some meta data for plots
ds_hydro_org = ds_hydro_org.raster.gdal_compliant() # spatial
ds_hydro_org["elevtn"].attrs.update(units="m", long_name="Elevation")
# plot elevation
# create nice elevation colormap
cmap_dem = colors.LinearSegmentedColormap.from_list("dem", plt.cm.terrain(np.linspace(0.25, 1, 256)))
norm_dem = colors.Normalize(vmin=0, vmax=2000)
ds_hydro_org["elevtn"].plot(
ax=ax,
zorder=1,
cbar_kwargs=dict(aspect=30, shrink=0.5),
alpha=0.5,
cmap=cmap_dem,
norm=norm_dem
)
# plot river network from original flow directions
gdf_riv_org.to_crs(da_flwdir.raster.crs).plot(
ax=ax, color="blue", linewidth=gdf_riv_org["strord"] / 3, label="Original flow directions"
)
# plot river network from new flow directions
gdf_riv.plot(
ax=ax, color="green", linewidth=gdf_riv["strord"] / 3, label="Derived flow directions",
)
ax.set_title("MERIT Hydro derived vs. original flow directions")
ax.legend()
[4]:
<matplotlib.legend.Legend at 0x7ff2e3dc6c90>
Deriving other DEM and flow directions related data#
Once you are satisfied with your flow direction map, you can create additionnal derived variables like upstream area or streamorder that can prove useful for example to build a model based on subbasin
region.
Note that to calculating upstream area and stream order requires all upstream areas to be included, which is not the case for all basins in this example. Here are some examples how to do that using PyFlwdir methods.
[5]:
# Create a new ds_hydro dataset with the riverburn flow directions
ds_hydro = da_flwdir.to_dataset(name="flwdir")
ds_hydro = ds_hydro.raster.gdal_compliant() # update spatial metadata
dims = ds_hydro.raster.dims
# add hydrological corrected elevation based on Yamazaki et al. (2012)
elevtn = flwdir.dem_adjust(elevtn=ds_hydro_org["elevtn"].values)
attrs = dict(_FillValue=-9999, long_name="corrected elevation", units="m")
ds_hydro["elevtn"] = xr.Variable(dims, elevtn, attrs=attrs)
# uparea (Note that this requires all upstream areas to be included.)
uparea = flwdir.upstream_area(unit="km2")
attrs = dict(_FillValue=-9999, long_name="upstream area", units="km2")
ds_hydro["uparea"] = xr.Variable(dims, uparea, attrs=attrs)
# stream order (Note that this requires all upstream areas to be included.)
strord = flwdir.stream_order()
attrs = dict(_FillValue=np.uint8(225), long_name="stream order", units="-")
ds_hydro["strord"] = xr.Variable(dims, strord, attrs=attrs)
# slope
slope = pyflwdir.dem.slope(
elevtn=ds_hydro["elevtn"].values,
nodata=ds_hydro["elevtn"].raster.nodata,
latlon=ds_hydro.raster.crs.is_geographic, # True if geographic crs, False if projected crs
transform=ds_hydro["elevtn"].raster.transform,
)
attrs = dict(_FillValue=-9999, long_name="slope", units="m/m")
ds_hydro["slope"] = xr.Variable(dims, slope, attrs=attrs)
# basin at the pits locations
basins = flwdir.basins(idxs=flwdir.idxs_pit).astype(np.int32)
attrs = dict(_FillValue=0, long_name="basin ids", units="-")
ds_hydro["basins"] = xr.Variable(dims, basins, attrs=attrs)
# basin index file
gdf_basins = ds_hydro["basins"].raster.vectorize()
ds_hydro
[5]:
<xarray.Dataset> Size: 2MB Dimensions: (latitude: 240, longitude: 360) Coordinates: * latitude (latitude) float64 2kB 46.2 46.2 46.2 46.2 ... 46.0 46.0 46.0 * longitude (longitude) float64 3kB 12.0 12.0 12.0 12.0 ... 12.3 12.3 12.3 spatial_ref int64 8B 0 Data variables: flwdir (latitude, longitude) uint8 86kB 0 16 1 0 16 1 ... 16 1 0 1 0 elevtn (latitude, longitude) float32 346kB 1.021e+03 ... 157.6 uparea (latitude, longitude) float64 691kB 0.2199 0.005943 ... 0.2386 strord (latitude, longitude) uint8 86kB 2 1 1 2 1 1 2 ... 2 1 1 3 1 2 slope (latitude, longitude) float32 346kB 0.5143 0.5413 ... 0.1146 basins (latitude, longitude) int32 346kB 1 1 2 2 2 ... 70 71 71 72 72
- latitude: 240
- longitude: 360
- latitude(latitude)float6446.2 46.2 46.2 ... 46.0 46.0 46.0
- standard_name :
- latitude
- long_name :
- latitude coordinate
- units :
- degrees_north
- axis :
- Y
array([46.199583, 46.19875 , 46.197917, ..., 46.002083, 46.00125 , 46.000417])
- longitude(longitude)float6412.0 12.0 12.0 ... 12.3 12.3 12.3
- standard_name :
- longitude
- long_name :
- longitude coordinate
- units :
- degrees_east
- axis :
- X
array([12.000417, 12.00125 , 12.002083, ..., 12.297917, 12.29875 , 12.299583])
- spatial_ref()int640
- crs_wkt :
- GEOGCRS["WGS 84",DATUM["World Geodetic System 1984",ELLIPSOID["WGS 84",6378137,298.257223563,LENGTHUNIT["metre",1]]],PRIMEM["Greenwich",0,ANGLEUNIT["degree",0.0174532925199433]],CS[ellipsoidal,2],AXIS["geodetic latitude (Lat)",north,ORDER[1],ANGLEUNIT["degree",0.0174532925199433]],AXIS["geodetic longitude (Lon)",east,ORDER[2],ANGLEUNIT["degree",0.0174532925199433]],ID["EPSG",4326]]
- semi_major_axis :
- 6378137.0
- semi_minor_axis :
- 6356752.314245179
- inverse_flattening :
- 298.257223563
- reference_ellipsoid_name :
- WGS 84
- longitude_of_prime_meridian :
- 0.0
- prime_meridian_name :
- Greenwich
- geographic_crs_name :
- WGS 84
- horizontal_datum_name :
- World Geodetic System 1984
- grid_mapping_name :
- latitude_longitude
- spatial_ref :
- GEOGCRS["WGS 84",DATUM["World Geodetic System 1984",ELLIPSOID["WGS 84",6378137,298.257223563,LENGTHUNIT["metre",1]]],PRIMEM["Greenwich",0,ANGLEUNIT["degree",0.0174532925199433]],CS[ellipsoidal,2],AXIS["geodetic latitude (Lat)",north,ORDER[1],ANGLEUNIT["degree",0.0174532925199433]],AXIS["geodetic longitude (Lon)",east,ORDER[2],ANGLEUNIT["degree",0.0174532925199433]],ID["EPSG",4326]]
- x_dim :
- longitude
- y_dim :
- latitude
- GeoTransform :
- 12.000000000000016 0.0008333333333333225 0.0 46.199999999999996 0.0 -0.0008333333333333186
array(0)
- flwdir(latitude, longitude)uint80 16 1 0 16 1 1 ... 16 16 1 0 1 0
- _FillValue :
- 247
- grid_mapping :
- spatial_ref
array([[ 0, 16, 1, ..., 2, 2, 4], [ 64, 32, 128, ..., 2, 2, 4], [ 64, 32, 32, ..., 4, 8, 8], ..., [ 64, 32, 128, ..., 8, 2, 4], [ 64, 32, 32, ..., 4, 2, 4], [ 64, 32, 32, ..., 0, 1, 0]], dtype=uint8)
- elevtn(latitude, longitude)float321.021e+03 1.073e+03 ... 157.6 157.6
- _FillValue :
- -9999
- long_name :
- corrected elevation
- units :
- m
array([[1020.8 , 1073.1 , 1105. , ..., 386.2 , 388.80002, 386.7 ], [1088.9 , 1122.8 , 1160.8 , ..., 386.2 , 384.7 , 383.80002], [1148.9 , 1171.8 , 1208.8 , ..., 386.2 , 381.7 , 380.1 ], ..., [ 793.8 , 811.10004, 824. , ..., 194.7 , 200.5 , 189.40001], [ 804.9 , 829.60004, 846.5 , ..., 189.90001, 189.5 , 177.40001], [ 823.5 , 845.3 , 859.5 , ..., 157.6 , 157.6 , 157.6 ]], dtype=float32)
- uparea(latitude, longitude)float640.2199 0.005943 ... 0.005965 0.2386
- _FillValue :
- -9999
- long_name :
- upstream area
- units :
- km2
array([[0.21991321, 0.00594301, 0.00594301, ..., 0.00594301, 0.00594301, 0.00594301], [0.19614088, 0.0118863 , 0.0059431 , ..., 0.0059431 , 0.01188612, 0.01782913], [0.1783113 , 0.01188648, 0.00594319, ..., 0.00594319, 0.0118863 , 0.03565844], ..., [0.03578676, 0.01789329, 0.00596434, ..., 0.00596434, 0.00596434, 0.20873662], [0.01789347, 0.01192895, 0.01192895, ..., 0.00596443, 0.00596443, 0.22066539], [0.00596452, 0.00596452, 0.00596452, ..., 0.10735596, 0.00596452, 0.23855886]])
- strord(latitude, longitude)uint82 1 1 2 1 1 2 3 ... 3 3 2 1 1 3 1 2
- _FillValue :
- 225
- long_name :
- stream order
- units :
- -
array([[2, 1, 1, ..., 1, 1, 1], [2, 1, 1, ..., 1, 1, 2], [2, 1, 1, ..., 1, 1, 2], ..., [2, 1, 1, ..., 1, 1, 2], [2, 1, 1, ..., 1, 1, 2], [1, 1, 1, ..., 3, 1, 2]], dtype=uint8)
- slope(latitude, longitude)float320.5143 0.5413 ... 0.1583 0.1146
- _FillValue :
- -9999
- long_name :
- slope
- units :
- m/m
array([[0.51432776, 0.5413488 , 0.27221325, ..., 0.00747009, 0.02149446, 0.01136101], [0.5459833 , 0.80539113, 0.6450317 , ..., 0.02068615, 0.03458828, 0.02887675], [0.474031 , 0.8115524 , 0.6958318 , ..., 0.05807183, 0.01134925, 0.02159704], ..., [0.17576262, 0.2855775 , 0.30759504, ..., 0.12126302, 0.11745112, 0.11047088], [0.22475238, 0.34088627, 0.27726272, ..., 0.17765708, 0.21685791, 0.15327132], [0.10499322, 0.23972316, 0.17005979, ..., 0.18411785, 0.1582703 , 0.1145691 ]], dtype=float32)
- basins(latitude, longitude)int321 1 2 2 2 3 3 ... 70 70 71 71 72 72
- _FillValue :
- 0
- long_name :
- basin ids
- units :
- -
array([[ 1, 1, 2, ..., 21, 21, 21], [ 1, 1, 2, ..., 21, 21, 21], [ 1, 1, 1, ..., 21, 21, 21], ..., [48, 48, 48, ..., 71, 72, 72], [48, 48, 48, ..., 71, 72, 72], [48, 48, 48, ..., 71, 72, 72]], dtype=int32)
- latitudePandasIndex
PandasIndex(Index([ 46.19958333333333, 46.19875, 46.197916666666664, 46.19708333333333, 46.19625, 46.19541666666667, 46.19458333333333, 46.193749999999994, 46.19291666666666, 46.19208333333333, ... 46.00791666666667, 46.00708333333333, 46.006249999999994, 46.00541666666666, 46.00458333333333, 46.00375, 46.002916666666664, 46.00208333333333, 46.00125, 46.000416666666666], dtype='float64', name='latitude', length=240))
- longitudePandasIndex
PandasIndex(Index([12.000416666666682, 12.001250000000017, 12.00208333333335, 12.002916666666682, 12.003750000000016, 12.004583333333349, 12.005416666666683, 12.006250000000016, 12.007083333333348, 12.007916666666683, ... 12.292083333333345, 12.292916666666677, 12.293750000000012, 12.294583333333344, 12.295416666666679, 12.296250000000011, 12.297083333333344, 12.297916666666678, 12.29875000000001, 12.299583333333345], dtype='float64', name='longitude', length=360))
[6]:
# NOTE that *upstream area* and *stream order* are not correct as not all upstream area is fully included in the domain.
# plot
fig, axes = plt.subplots(2, 2, figsize=(12, 10), sharex=True, sharey=True)
# plot uparea; use a log scale colormap
ds_hydro["uparea"].plot(
ax=axes[0, 0],
norm=colors.LogNorm(vmin=1, vmax=ds_hydro["uparea"].max()),
cmap='viridis',
)
axes[0,0].set_title("Upstream area [km2]")
# plot strord
ds_hydro["strord"].plot(
ax=axes[0, 1],
cmap=colors.ListedColormap(cm.Blues(np.linspace(0, 1, 7)))
)
axes[0,1].set_title("Strahler Stream order")
# plot slope
ds_hydro["slope"].plot(ax=axes[1, 1])
axes[1,1].set_title("Slope")
# plot basins
# elevation cmap
ds_hydro["elevtn"].plot(ax=axes[1, 0], cmap=cmap_dem, norm=norm_dem)
axes[1,0].set_title("Elevation")
[6]:
Text(0.5, 1.0, 'Elevation')
Exporting the newly created data and corresponding data catalog#
Finally, once we are happy with the new dataset, we can write out the data and create the corresponding data catalog so that it can be re-used to build a new wflow model.
[7]:
# Export the gridded data as tif files in a new folder
output_path = "./hydro_data"
# export the hydrography data as tif files (one per variable)
ds_hydro.raster.to_mapstack(
root=os.path.join(output_path, "ds_hydro"),
driver="GTiff",
)
# export the basin index as geosjon
gdf_basins.to_file(
os.path.join(output_path, "da_hydro_basins.geojson"), driver="GeoJSON"
)
Now let’s prepare the corresponding data catalog: (the writefile command will directly write a file using the lines in the jupyter cell)
[8]:
%%writefile ./hydro_data/data_catalog.yml
ds_hydro:
data_type: RasterDataset
driver: raster
path: ./ds_hydro/{variable}.tif
rename:
slope: lndslp
meta:
category: topography
processing_notes: prepared from MERIT Hydro using hydromt d8_from_dem and pyflwdir
processing_script: prepare_ldd.ipynb from hydromt_wflow repository
da_hydro_new_index:
data_type: GeoDataFrame
driver: vector
path: ./da_hydro_basins.geojson
meta:
processing_notes: prepared from MERIT Hydro using hydromt d8_from_dem and pyflwdir
processing_script: prepare_ldd.ipynb from hydromt_wflow repository
Writing ./hydro_data/data_catalog.yml
And now let’s try to load our data again with hydromt:
[9]:
data_catalog = DataCatalog(data_libs="./hydro_data/data_catalog.yml")
ds_hydro_new = data_catalog.get_rasterdataset("ds_hydro")
ds_hydro_new
[9]:
<xarray.Dataset> Size: 2MB Dimensions: (x: 360, y: 240) Coordinates: * x (x) float64 3kB 12.0 12.0 12.0 12.0 ... 12.3 12.3 12.3 12.3 * y (y) float64 2kB 46.2 46.2 46.2 46.2 ... 46.0 46.0 46.0 46.0 spatial_ref int64 8B 0 Data variables: basins (y, x) int32 346kB dask.array<chunksize=(240, 360), meta=np.ndarray> elevtn (y, x) float32 346kB dask.array<chunksize=(240, 360), meta=np.ndarray> flwdir (y, x) uint8 86kB dask.array<chunksize=(240, 360), meta=np.ndarray> lndslp (y, x) float32 346kB dask.array<chunksize=(240, 360), meta=np.ndarray> strord (y, x) uint8 86kB dask.array<chunksize=(240, 360), meta=np.ndarray> uparea (y, x) float64 691kB dask.array<chunksize=(240, 360), meta=np.ndarray> Attributes: category: topography processing_notes: prepared from MERIT Hydro using hydromt d8_from_dem a... processing_script: prepare_ldd.ipynb from hydromt_wflow repository
- x: 360
- y: 240
- x(x)float6412.0 12.0 12.0 ... 12.3 12.3 12.3
array([12.000417, 12.00125 , 12.002083, ..., 12.297917, 12.29875 , 12.299583])
- y(y)float6446.2 46.2 46.2 ... 46.0 46.0 46.0
array([46.199583, 46.19875 , 46.197917, ..., 46.002083, 46.00125 , 46.000417])
- spatial_ref()int640
- crs_wkt :
- GEOGCRS["WGS 84",DATUM["World Geodetic System 1984",ELLIPSOID["WGS 84",6378137,298.257223563,LENGTHUNIT["metre",1]]],PRIMEM["Greenwich",0,ANGLEUNIT["degree",0.0174532925199433]],CS[ellipsoidal,2],AXIS["geodetic latitude (Lat)",north,ORDER[1],ANGLEUNIT["degree",0.0174532925199433]],AXIS["geodetic longitude (Lon)",east,ORDER[2],ANGLEUNIT["degree",0.0174532925199433]],ID["EPSG",4326]]
- semi_major_axis :
- 6378137.0
- semi_minor_axis :
- 6356752.314245179
- inverse_flattening :
- 298.257223563
- reference_ellipsoid_name :
- WGS 84
- longitude_of_prime_meridian :
- 0.0
- prime_meridian_name :
- Greenwich
- geographic_crs_name :
- WGS 84
- horizontal_datum_name :
- World Geodetic System 1984
- grid_mapping_name :
- latitude_longitude
- spatial_ref :
- GEOGCRS["WGS 84",DATUM["World Geodetic System 1984",ELLIPSOID["WGS 84",6378137,298.257223563,LENGTHUNIT["metre",1]]],PRIMEM["Greenwich",0,ANGLEUNIT["degree",0.0174532925199433]],CS[ellipsoidal,2],AXIS["geodetic latitude (Lat)",north,ORDER[1],ANGLEUNIT["degree",0.0174532925199433]],AXIS["geodetic longitude (Lon)",east,ORDER[2],ANGLEUNIT["degree",0.0174532925199433]],ID["EPSG",4326]]
- GeoTransform :
- 12.000000000000016 0.0008333333333333225 0.0 46.199999999999996 0.0 -0.0008333333333333186
- x_dim :
- x
- y_dim :
- y
array(0)
- basins(y, x)int32dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- 0
- source_file :
- basins.tif
Array Chunk Bytes 337.50 kiB 337.50 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type int32 numpy.ndarray - elevtn(y, x)float32dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- -9999.0
- source_file :
- elevtn.tif
Array Chunk Bytes 337.50 kiB 337.50 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type float32 numpy.ndarray - flwdir(y, x)uint8dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- 247
- source_file :
- flwdir.tif
Array Chunk Bytes 84.38 kiB 84.38 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type uint8 numpy.ndarray - lndslp(y, x)float32dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- -9999.0
- source_file :
- slope.tif
Array Chunk Bytes 337.50 kiB 337.50 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type float32 numpy.ndarray - strord(y, x)uint8dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- 225
- source_file :
- strord.tif
Array Chunk Bytes 84.38 kiB 84.38 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type uint8 numpy.ndarray - uparea(y, x)float64dask.array<chunksize=(240, 360), meta=np.ndarray>
- AREA_OR_POINT :
- Area
- _FillValue :
- -9999.0
- source_file :
- uparea.tif
Array Chunk Bytes 675.00 kiB 675.00 kiB Shape (240, 360) (240, 360) Dask graph 1 chunks in 3 graph layers Data type float64 numpy.ndarray
- xPandasIndex
PandasIndex(Index([12.000416666666682, 12.001250000000015, 12.00208333333335, 12.002916666666682, 12.003750000000016, 12.004583333333349, 12.005416666666683, 12.006250000000016, 12.007083333333348, 12.007916666666683, ... 12.292083333333345, 12.29291666666668, 12.293750000000012, 12.294583333333344, 12.295416666666679, 12.296250000000011, 12.297083333333346, 12.297916666666678, 12.298750000000013, 12.299583333333345], dtype='float64', name='x', length=360))
- yPandasIndex
PandasIndex(Index([ 46.19958333333333, 46.19875, 46.197916666666664, 46.19708333333333, 46.19625, 46.19541666666666, 46.19458333333333, 46.193749999999994, 46.19291666666666, 46.19208333333333, ... 46.00791666666667, 46.007083333333334, 46.00625, 46.00541666666667, 46.004583333333336, 46.00375, 46.002916666666664, 46.00208333333333, 46.00125, 46.000416666666666], dtype='float64', name='y', length=240))
- category :
- topography
- processing_notes :
- prepared from MERIT Hydro using hydromt d8_from_dem and pyflwdir
- processing_script :
- prepare_ldd.ipynb from hydromt_wflow repository