xarray ODIM backend#

In this example, we read ODIM_H5 (HDF5) data files using the xradar odim backend. Throughout the notebook xarray accessors are used to access wradlib functionality.

[1]:
import glob
import os
import wradlib as wrl
import warnings

warnings.filterwarnings("ignore")
import matplotlib.pyplot as plt
import numpy as np
import xradar as xd
import datatree as xt
import xarray as xr

try:
    get_ipython().run_line_magic("matplotlib inline")
except:
    plt.ion()
/home/runner/micromamba/envs/wradlib-tests/lib/python3.11/site-packages/h5py/__init__.py:36: UserWarning: h5py is running against HDF5 1.14.3 when it was built against 1.14.2, this may cause problems
  _warn(("h5py is running against HDF5 {0} when it was built against {1}, "

Load ODIM_H5 Volume Data#

[2]:
fpath = "hdf5/knmi_polar_volume.h5"
f = wrl.util.get_wradlib_data_file(fpath)
vol = xd.io.open_odim_datatree(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)
<xarray.DatasetView>
Dimensions:              ()
Data variables:
    volume_number        int64 0
    platform_type        <U5 'fixed'
    instrument_type      <U5 'radar'
    time_coverage_start  <U20 '2011-06-10T11:40:02Z'
    time_coverage_end    <U20 '2011-06-10T11:43:54Z'
    longitude            float32 4.79
    altitude             float32 50.0
    latitude             float32 52.95
Attributes:
    Conventions:      ODIM_H5/V2_0
    version:          None
    title:            None
    institution:      None
    references:       None
    source:           None
    history:          None
    comment:          im/exported using xradar
    instrument_name:  None

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.DatasetView>
Dimensions:              ()
Data variables:
    volume_number        int64 0
    platform_type        <U5 'fixed'
    instrument_type      <U5 'radar'
    time_coverage_start  <U20 '2011-06-10T11:40:02Z'
    time_coverage_end    <U20 '2011-06-10T11:43:54Z'
    longitude            float32 4.79
    altitude             float32 50.0
    latitude             float32 52.95
Attributes:
    Conventions:      ODIM_H5/V2_0
    version:          None
    title:            None
    institution:      None
    references:       None
    source:           None
    history:          None
    comment:          im/exported using xradar
    instrument_name:  None

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["sweep_0"])
<xarray.DatasetView>
Dimensions:            (azimuth: 360, range: 320)
Coordinates:
    elevation          (azimuth) float32 ...
    time               (azimuth) datetime64[ns] 2011-06-10T11:40:17.361118208...
  * range              (range) float32 500.0 1.5e+03 ... 3.185e+05 3.195e+05
    longitude          float32 ...
    latitude           float32 ...
    altitude           float32 ...
  * azimuth            (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5
Data variables:
    DBZH               (azimuth, range) float32 ...
    sweep_mode         <U20 ...
    sweep_number       int64 ...
    prt_mode           <U7 ...
    follow_mode        <U7 ...
    sweep_fixed_angle  float32 ...

Georeferencing#

[6]:
swp = vol["sweep_0"].ds
swp = swp.assign_coords(sweep_mode=swp.sweep_mode)
swp = swp.wrl.georef.georeference()

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.

[7]:
display(swp.DBZH)
<xarray.DataArray 'DBZH' (azimuth: 360, range: 320)>
[115200 values with dtype=float32]
Coordinates: (12/15)
    elevation   (azimuth) float32 0.3 0.3 0.3 0.3 0.3 ... 0.3 0.3 0.3 0.3 0.3
    time        (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
    latitude    float32 52.95
    ...          ...
    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
    crs_wkt     int64 0
Attributes:
    _Undetect:      0.0
    units:          dBZ
    standard_name:  radar_equivalent_reflectivity_factor_h
    long_name:      Equivalent reflectivity factor H

Create simple plot#

Using xarray features a simple plot can be created like this. Note the sortby('time') method, which sorts the radials by time.

For more details on plotting radar data see under Visualization.

[8]:
swp.DBZH.sortby("time").plot(x="range", y="time", add_labels=False)
[8]:
<matplotlib.collections.QuadMesh at 0x7f0b4ff5bc50>
../../../_images/notebooks_fileio_backends_odim_backend_16_1.png
[9]:
fig = plt.figure(figsize=(10, 10))
pm = swp.DBZH.wrl.vis.plot(crs={"latmin": 33e3}, fig=fig)
../../../_images/notebooks_fileio_backends_odim_backend_17_0.png

Retrieve explicit group#

[10]:
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:
    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 ...
  * azimuth            (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5
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#

[11]:
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="odim", concat_dim="volume_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:            (volume_time: 2, azimuth: 360, range: 1200)
Coordinates:
    elevation          (azimuth) float32 dask.array<chunksize=(360,), meta=np.ndarray>
    time               (volume_time, azimuth) datetime64[ns] 2018-12-20T06:06...
  * range              (range) float32 125.0 375.0 625.0 ... 2.996e+05 2.999e+05
    longitude          float64 151.2
    latitude           float64 -33.7
    altitude           float64 195.0
  * azimuth            (azimuth) float32 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5
Dimensions without coordinates: volume_time
Data variables: (12/17)
    DBZH               (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    DBZH_CLEAN         (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    VRADDH             (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    VRADH              (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    WRADH              (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    TH                 (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    ...                 ...
    CLASS              (volume_time, azimuth, range) float32 dask.array<chunksize=(1, 360, 1200), meta=np.ndarray>
    sweep_mode         (volume_time) <U20 'azimuth_surveillance' 'azimuth_sur...
    sweep_number       (volume_time) int64 0 0
    prt_mode           (volume_time) <U7 'not_set' 'not_set'
    follow_mode        (volume_time) <U7 'not_set' 'not_set'
    sweep_fixed_angle  (volume_time) float64 0.5 0.5
[ ]: