xarray ODIM backend#
In this example, we read ODIM_H5 (HDF5) data files using the xarray odim
backend.
[1]:
import glob
import os
import wradlib as wrl
import warnings
warnings.filterwarnings("ignore")
import matplotlib.pyplot as pl
import numpy as np
import xarray as xr
try:
get_ipython().run_line_magic("matplotlib inline")
except:
pl.ion()
from wradlib.io import open_odim_dataset
/home/runner/micromamba-root/envs/wradlib-tests/lib/python3.11/site-packages/tqdm/auto.py:22: TqdmWarning: IProgress not found. Please update jupyter and ipywidgets. See https://ipywidgets.readthedocs.io/en/stable/user_install.html
from .autonotebook import tqdm as notebook_tqdm
Load ODIM_H5 Volume Data#
[2]:
fpath = "hdf5/knmi_polar_volume.h5"
f = wrl.util.get_wradlib_data_file(fpath)
vol = wrl.io.open_odim_dataset(f)
Downloading file 'hdf5/knmi_polar_volume.h5' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/hdf5/knmi_polar_volume.h5' to '/home/runner/work/wradlib/wradlib/wradlib-data'.
Inspect RadarVolume#
[3]:
display(vol)
<wradlib.RadarVolume>
Dimension(s): (sweep: 14)
Elevation(s): (0.3, 0.4, 0.8, 1.1, 2.0, 3.0, 4.5, 6.0, 8.0, 10.0, 12.0, 15.0, 20.0, 25.0)
Inspect root group#
The sweep
dimension contains the number of scans in this radar volume. Further the dataset consists of variables (location coordinates, time_coverage) and attributes (Conventions, metadata).
[4]:
vol.root
[4]:
<xarray.Dataset> Dimensions: (sweep: 14) Dimensions without coordinates: sweep Data variables: volume_number int64 0 platform_type <U5 'fixed' instrument_type <U5 'radar' primary_axis <U6 'axis_z' time_coverage_start <U20 '2011-06-10T11:40:02Z' time_coverage_end <U20 '2011-06-10 11:43:55Z' latitude float32 ... longitude float32 ... altitude float32 ... sweep_group_name (sweep) <U8 'sweep_0' 'sweep_1' ... 'sweep_13' sweep_fixed_angle (sweep) float32 0.3 0.4 0.8 1.1 ... 12.0 15.0 20.0 25.0 Attributes: version: None title: None institution: None references: None source: None history: None comment: im/exported using wradlib instrument_name: None fixed_angle: 0.3
Inspect sweep group(s)#
The sweep-groups can be accessed via their respective keys. The dimensions consist of range
and time
with added coordinates azimuth
, elevation
, range
and time
. There will be variables like radar moments (DBZH etc.) and sweep-dependend metadata (like fixed_angle
, sweep_mode
etc.).
[5]:
display(vol[0])
<xarray.Dataset> Dimensions: (azimuth: 360, range: 320) Coordinates: * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5 elevation (azimuth) float32 ... rtime (azimuth) datetime64[ns] 2011-06-10T11:40:17.361118208... * range (range) float32 500.0 1.5e+03 ... 3.185e+05 3.195e+05 sweep_mode <U20 ... longitude float32 ... latitude float32 ... altitude float32 ... time datetime64[ns] 2011-06-10T11:40:02.027777792 Data variables: DBZH (azimuth, range) float32 ... sweep_number int64 ... prt_mode <U7 ... follow_mode <U7 ... sweep_fixed_angle float32 0.3 Attributes: fixed_angle: 0.3
Goereferencing#
[6]:
swp = vol[0].copy().pipe(wrl.georef.georeference_dataset)
Plotting#
[7]:
swp.DBZH.plot.pcolormesh(x="x", y="y")
pl.gca().set_aspect("equal")

[8]:
fig = pl.figure(figsize=(10, 10))
swp.DBZH.wradlib.plot_ppi(proj="cg", fig=fig)
[8]:
<matplotlib.collections.QuadMesh at 0x7fa79317b0d0>

[9]:
import cartopy
import cartopy.crs as ccrs
import cartopy.feature as cfeature
map_trans = ccrs.AzimuthalEquidistant(
central_latitude=swp.latitude.values, central_longitude=swp.longitude.values
)
[10]:
map_proj = ccrs.AzimuthalEquidistant(
central_latitude=swp.latitude.values, central_longitude=swp.longitude.values
)
pm = swp.DBZH.wradlib.plot_ppi(proj=map_proj)
ax = pl.gca()
ax.gridlines(crs=map_proj)
print(ax)
< GeoAxes: +proj=aeqd +ellps=WGS84 +lon_0=4.7899699211120605 +lat_0=52.953338623046875 +x_0=0.0 +y_0=0.0 +no_defs +type=crs >

[11]:
map_proj = ccrs.Mercator(central_longitude=swp.longitude.values)
fig = pl.figure(figsize=(10, 8))
ax = fig.add_subplot(111, projection=map_proj)
pm = swp.DBZH.wradlib.plot_ppi(ax=ax)
ax.gridlines(draw_labels=True)
[11]:
<cartopy.mpl.gridliner.Gridliner at 0x7fa78ac4ab10>

[12]:
import cartopy.feature as cfeature
def plot_borders(ax):
borders = cfeature.NaturalEarthFeature(
category="physical", name="coastline", scale="10m", facecolor="none"
)
ax.add_feature(borders, edgecolor="black", lw=2, zorder=4)
map_proj = ccrs.Mercator(central_longitude=swp.longitude.values)
fig = pl.figure(figsize=(10, 8))
ax = fig.add_subplot(111, projection=map_proj)
DBZH = swp.DBZH
pm = DBZH.where(DBZH > 0).wradlib.plot_ppi(ax=ax)
plot_borders(ax)
ax.gridlines(draw_labels=True)
[12]:
<cartopy.mpl.gridliner.Gridliner at 0x7fa78adef290>

[13]:
import matplotlib.path as mpath
theta = np.linspace(0, 2 * np.pi, 100)
center, radius = [0.5, 0.5], 0.5
verts = np.vstack([np.sin(theta), np.cos(theta)]).T
circle = mpath.Path(verts * radius + center)
map_proj = ccrs.AzimuthalEquidistant(
central_latitude=swp.latitude.values,
central_longitude=swp.longitude.values,
)
fig = pl.figure(figsize=(10, 8))
ax = fig.add_subplot(111, projection=map_proj)
ax.set_boundary(circle, transform=ax.transAxes)
pm = swp.DBZH.wradlib.plot_ppi(proj=map_proj, ax=ax)
ax = pl.gca()
ax.gridlines(crs=map_proj)
[13]:
<cartopy.mpl.gridliner.Gridliner at 0x7fa78b025610>

[14]:
fig = pl.figure(figsize=(10, 8))
proj = ccrs.AzimuthalEquidistant(
central_latitude=swp.latitude.values, central_longitude=swp.longitude.values
)
ax = fig.add_subplot(111, projection=proj)
pm = swp.DBZH.wradlib.plot_ppi(ax=ax)
ax.gridlines()
[14]:
<cartopy.mpl.gridliner.Gridliner at 0x7fa78afa3bd0>

[15]:
swp.DBZH.wradlib.plot_ppi()
[15]:
<matplotlib.collections.QuadMesh at 0x7fa7886a74d0>

Inspect radar moments#
The DataArrays can be accessed by key or by attribute. Each DataArray has dimensions and coordinates of it’s parent dataset. There are attributes connected which are defined by ODIM_H5 standard.
[16]:
display(swp.DBZH)
<xarray.DataArray 'DBZH' (azimuth: 360, range: 320)> array([[ 22. , 17. , -8. , ..., -31.5, -31.5, -31.5], [ 24. , 24.5, -9. , ..., -31.5, -31.5, -31.5], [ 35.5, 42. , 12. , ..., -31.5, -31.5, -31.5], ..., [ 23. , 14. , -13. , ..., -31.5, -31.5, -31.5], [ 23. , 14. , -9. , ..., -31.5, -31.5, -31.5], [ 22. , 18.5, -11.5, ..., -31.5, -31.5, -31.5]], dtype=float32) Coordinates: (12/15) * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 356.5 357.5 358.5 359.5 elevation (azimuth) float32 0.3 0.3 0.3 0.3 0.3 ... 0.3 0.3 0.3 0.3 0.3 rtime (azimuth) datetime64[ns] 2011-06-10T11:40:17.361118208 ... 20... * range (range) float32 500.0 1.5e+03 2.5e+03 ... 3.185e+05 3.195e+05 sweep_mode <U20 'azimuth_surveillance' longitude float32 4.79 ... ... x (azimuth, range) float32 4.363 13.09 ... -2.777e+03 -2.786e+03 y (azimuth, range) float32 500.0 1.5e+03 ... 3.183e+05 3.193e+05 z (azimuth, range) float32 53.0 58.0 64.0 ... 7.691e+03 7.734e+03 gr (azimuth, range) float32 500.0 1.5e+03 ... 3.183e+05 3.193e+05 rays (azimuth, range) float32 0.5 0.5 0.5 0.5 ... 359.5 359.5 359.5 bins (azimuth, range) float32 500.0 1.5e+03 ... 3.185e+05 3.195e+05 Attributes: _Undetect: 0.0 long_name: Equivalent reflectivity factor H units: dBZ standard_name: radar_equivalent_reflectivity_factor_h
Create simple plot#
Using xarray features a simple plot can be created like this. Note the sortby('rtime')
method, which sorts the radials by time.
[17]:
swp.DBZH.sortby("rtime").plot(x="range", y="rtime", add_labels=False)
[17]:
<matplotlib.collections.QuadMesh at 0x7fa793112050>

[18]:
fig = pl.figure(figsize=(5, 5))
pm = swp.DBZH.wradlib.plot_ppi(proj={"latmin": 33e3}, fig=fig)

Mask some values#
[19]:
swp["DBZH"] = swp["DBZH"].where(swp["DBZH"] >= 0)
swp["DBZH"].plot()
[19]:
<matplotlib.collections.QuadMesh at 0x7fa788555690>

Export to ODIM and CfRadial2#
[20]:
vol.to_odim("knmi_odim.h5")
[21]:
vol.to_cfradial2("knmi_odim_as_cfradial2.nc")
Import again#
[22]:
vola = wrl.io.open_odim_dataset("knmi_odim.h5")
[23]:
volb = wrl.io.open_cfradial2_dataset("knmi_odim_as_cfradial2.nc")
Check equality#
[24]:
xr.testing.assert_allclose(vol.root, vola.root)
xr.testing.assert_equal(vol[0], vola[0])
xr.testing.assert_allclose(vol.root, volb.root)
xr.testing.assert_equal(vol[0], volb[0])
xr.testing.assert_allclose(vola.root, volb.root)
xr.testing.assert_equal(vola[0], volb[0])
More ODIM loading mechanisms#
Use xr.open_dataset
to retrieve explicit group#
Warning
Since \(\omega radlib\) version 1.18 the xarray backend engines for polar radar data have been renamed and prepended with wradlib-
(eg. odim
-> wradlib-odim
). This was necessary to avoid clashes with the new xradar-package, which will eventually replace the wradlib engines. Users have to make sure to check which engine to use for their use-case when using xarray.open_dataset
. Users might install and test xradar
, and check if it is
already robust enough for their use-cases (by using xradar’s engine="odim"
.
Since \(\omega radlib\) version 1.19 the xarray backend engines for polar radar data have been deprecated. The functionality is kept until wradlib version 2.0, when the backend-code will be removed completely. wradlib is importing that functionality from xradar-package whenever and wherever necessary.
Below we use a compatibility layer in wradlib to give users the chance to adapt their code. The first minimal change is that for every backend the group-layout is conforming to the CfRadial-standard naming scheme (sweep_0
, sweep_1
, etc.).
Below you can inspect the main differences of the wradlib compatibility layer and the plain xradar implementation.
use wradlib compatibility layer#
[25]:
swp_a = xr.open_dataset(
f, engine="wradlib-odim", group="sweep_13", backend_kwargs=dict(reindex_angle=False)
)
display(swp_a)
<xarray.Dataset> Dimensions: (azimuth: 360, range: 240) Coordinates: * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5 elevation (azimuth) float32 ... rtime (azimuth) datetime64[ns] ... * range (range) float32 250.0 750.0 ... 1.192e+05 1.198e+05 sweep_mode <U20 ... longitude float32 ... latitude float32 ... altitude float32 ... time datetime64[ns] ... Data variables: DBZH (azimuth, range) float32 ... sweep_number int64 ... prt_mode <U7 ... follow_mode <U7 ... sweep_fixed_angle float32 ... Attributes: fixed_angle: 25.0
use xradar backend#
[26]:
swp_b = xr.open_dataset(
f, engine="odim", group="sweep_13", backend_kwargs=dict(reindex_angle=False)
)
display(swp_b)
<xarray.Dataset> Dimensions: (azimuth: 360, range: 240) Coordinates: * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5 elevation (azimuth) float32 ... time (azimuth) datetime64[ns] ... * range (range) float32 250.0 750.0 ... 1.192e+05 1.198e+05 longitude float32 ... latitude float32 ... altitude float32 ... Data variables: DBZH (azimuth, range) float32 ... sweep_mode <U20 ... sweep_number int64 ... prt_mode <U7 ... follow_mode <U7 ... sweep_fixed_angle float32 ...
Use xr.open_mfdataset
to retrieve timeseries of explicit group#
[27]:
flist = ["hdf5/71_20181220_060628.pvol.h5", "hdf5/71_20181220_061228.pvol.h5"]
flist = [wrl.util.get_wradlib_data_file(f) for f in flist]
ts = xr.open_mfdataset(
flist, engine="wradlib-odim", concat_dim="time", combine="nested", group="sweep_0"
)
display(ts)
Downloading file 'hdf5/71_20181220_060628.pvol.h5' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/hdf5/71_20181220_060628.pvol.h5' to '/home/runner/work/wradlib/wradlib/wradlib-data'.
Downloading file 'hdf5/71_20181220_061228.pvol.h5' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/hdf5/71_20181220_061228.pvol.h5' to '/home/runner/work/wradlib/wradlib/wradlib-data'.
<xarray.Dataset> Dimensions: (time: 2, azimuth: 360, range: 1200) Coordinates: * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5 elevation (azimuth) float32 dask.array<chunksize=(360,), meta=np.ndarray> rtime (time, azimuth) datetime64[ns] 2018-12-20T06:06:50.112... * range (range) float32 125.0 375.0 625.0 ... 2.996e+05 2.999e+05 sweep_mode <U20 'azimuth_surveillance' longitude float64 151.2 latitude float64 -33.7 altitude float64 195.0 * time (time) datetime64[ns] 2018-12-20T06:06:28.040277760 20... Data variables: (12/16) DBZH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> DBZH_CLEAN (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> VRADDH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> VRADH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> WRADH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> TH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> ... ... SNRH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> CLASS (time, azimuth, range) int8 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> sweep_number (time) int64 0 0 prt_mode (time) <U7 'not_set' 'not_set' follow_mode (time) <U7 'not_set' 'not_set' sweep_fixed_angle (time) float64 0.5 0.5 Attributes: fixed_angle: 0.5
Use wrl.io.open_odim_mfdataset
to retrieve volume timeseries#
[28]:
fpath = os.path.join(wrl.util.get_wradlib_data_path(), "hdf5/71*.h5")
f = glob.glob(fpath)
ts = wrl.io.open_odim_mfdataset(f)
display(ts)
100%|██████████| 14/14 [00:04<00:00, 3.24it/s]
[29]:
display(ts[0])
<xarray.Dataset> Dimensions: (time: 2, azimuth: 360, range: 1200) Coordinates: * azimuth (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5 elevation (azimuth) float32 dask.array<chunksize=(360,), meta=np.ndarray> rtime (time, azimuth) datetime64[ns] 2018-12-20T06:06:50.112... * range (range) float32 125.0 375.0 625.0 ... 2.996e+05 2.999e+05 sweep_mode <U20 'azimuth_surveillance' longitude float64 151.2 latitude float64 -33.7 altitude float64 195.0 * time (time) datetime64[ns] 2018-12-20T06:06:28.040277760 20... Data variables: (12/16) DBZH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> DBZH_CLEAN (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> VRADDH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> VRADH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> WRADH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> TH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> ... ... SNRH (time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> CLASS (time, azimuth, range) int8 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray> sweep_number (time) int64 0 0 prt_mode (time) <U7 'not_set' 'not_set' follow_mode (time) <U7 'not_set' 'not_set' sweep_fixed_angle (time) float64 0.5 0.5 Attributes: fixed_angle: 0.5
[ ]: