xarray powered Cf/Radial and ODIM_H5

In this example, we read and write Cf/Radial (NetCDF) and ODIM_H5 (HDF5) data files from different sources using an xarray powered data structure.

[1]:
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().magic("matplotlib inline")
except:
    pl.ion()
from wradlib.io.xarray import CfRadial, OdimH5

Load ODIM_H5 Volume Data

[2]:
fpath = 'hdf5/knmi_polar_volume.h5'
f = wrl.util.get_wradlib_data_file(fpath)
cf1 = OdimH5(f, standard='cf', georef=True)

Inspect root group

You can use the object dictionary using cf1[‘root’] or the property cf1.root.

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).

[3]:
cf1.root
[3]:
<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-10T11:43:54Z'
    latitude             float32 52.95334
    longitude            float32 4.78997
    altitude             float32 50.0
    altitude_agl         float64 nan
    sweep_group_name     (sweep) <U8 'sweep_1' 'sweep_2' ... 'sweep_14'
    sweep_fixed_angle    (sweep) float32 0.3 0.4 0.8 1.1 ... 12.0 15.0 20.0 25.0
    frequency            float64 nan
    status_xml           <U4 'None'
Attributes:
    Conventions:          Cf/Radial
    version:              H5rad 2.0
    title:                None
    institution:          RAD:NL51;PLC:nldhl
    references:           None
    source:               None
    history:              None
    comment:              im/exported using wradlib
    instrument_name:      None
    site_name:            name of site where data were gathered
    scan_name:            name of scan strategy used, if applicable
    scan_id:              scan strategy id, if applicable. assumed 0 if missing
    platform_is_mobile:   "true" or "false", assumed "false" if missing
    ray_times_increase:   "true" or "false", assumed "true" if missing. This ...
    field_names:          array of strings of field names present in this file.
    time_coverage_start:  copy of time_coverage_start global variable
    time_coverage_end:    copy of time_coverage_end global variable
    simulated data:       "true" or "false", assumed "false" if missing. data...
    instrument:           RAD:NL51;PLC:nldhl

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.).

[4]:
cf1['sweep_1']
[4]:
<xarray.Dataset>
Dimensions:       (range: 320, time: 360)
Coordinates:
    sweep_mode    <U20 ...
    latitude      float32 ...
    altitude      float32 ...
    longitude     float32 ...
    elevation     (time) float32 ...
    azimuth       (time) float32 ...
  * range         (range) float32 500.0 1500.0 2500.0 ... 318500.0 319500.0
    y             (time, range) float32 ...
    z             (time, range) float32 ...
    gr            (time, range) float32 ...
    rays          (time, range) float32 ...
    bins          (time, range) float32 ...
    x             (time, range) float32 ...
  * time          (time) datetime64[ns] 2011-06-10T11:40:06.694446592 ... 2011-06-10T11:40:06.638891008
Data variables:
    DBZH          (time, range) float32 ...
    sweep_number  int64 ...
    follow_mode   <U4 ...
    prt_mode      <U5 ...
    fixed_angle   float32 ...
[5]:
cf1['sweep_1'].DBZH
[5]:
<xarray.DataArray 'DBZH' (time: 360, range: 320)>
[115200 values with dtype=float32]
Coordinates:
    sweep_mode  <U20 ...
    latitude    float32 ...
    altitude    float32 ...
    longitude   float32 ...
    elevation   (time) float32 ...
    azimuth     (time) float32 ...
  * range       (range) float32 500.0 1500.0 2500.0 ... 318500.0 319500.0
    y           (time, range) float32 ...
    z           (time, range) float32 ...
    gr          (time, range) float32 ...
    rays        (time, range) float32 ...
    bins        (time, range) float32 ...
    x           (time, range) float32 ...
  * time        (time) datetime64[ns] 2011-06-10T11:40:06.694446592 ... 2011-06-10T11:40:06.638891008
Attributes:
    IMAGE_VERSION:  1.2
    standard_name:  radar_equivalent_reflectivity_factor_h
    long_name:      Equivalent reflectivity factor H
    units:          dBZ

Plotting

[6]:
cf1['sweep_1'].DBZH.plot.pcolormesh(x='x', y='y')
pl.gca().set_aspect('equal')
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_11_0.png
[7]:
fig = pl.figure(figsize=(10,8))
cf1['sweep_1'].DBZH.sortby('azimuth').wradlib.plot_ppi(proj='cg', fig=fig)
[7]:
<matplotlib.collections.QuadMesh at 0x7f49a38eefd0>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_12_1.png
[8]:
import cartopy
import cartopy.crs as ccrs
import cartopy.feature as cfeature

map_trans = ccrs.AzimuthalEquidistant(central_latitude=cf1['sweep_1'].latitude.values,
                                      central_longitude=cf1['sweep_1'].longitude.values)
[9]:
map_proj = ccrs.AzimuthalEquidistant(central_latitude=cf1['sweep_1'].latitude.values,
                                      central_longitude=cf1['sweep_1'].longitude.values)
pm = cf1['sweep_1'].DBZH.wradlib.plot_ppi(proj=map_proj)
ax = pl.gca()
ax.gridlines(crs=map_proj)
print(ax)
< GeoAxes: <cartopy.crs.AzimuthalEquidistant object at 0x7f49b16a7900> >
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_14_1.png
[10]:
map_proj = ccrs.Mercator(central_longitude=cf1['sweep_1'].longitude.values)
fig = pl.figure(figsize=(10,8))
ax = fig.add_subplot(111, projection=map_proj)
pm = cf1['sweep_1'].DBZH.wradlib.plot_ppi(ax=ax)
ax.gridlines(draw_labels=True)
[10]:
<cartopy.mpl.gridliner.Gridliner at 0x7f49b15c6f10>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_15_1.png
[11]:
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=cf1['sweep_1'].longitude.values)
fig = pl.figure(figsize=(10,8))
ax = fig.add_subplot(111, projection=map_proj)

DBZH = cf1['sweep_1'].DBZH
pm = DBZH.where(DBZH > 0).wradlib.plot_ppi(ax=ax)
plot_borders(ax)
ax.gridlines(draw_labels=True)
[11]:
<cartopy.mpl.gridliner.Gridliner at 0x7f49a0a90c40>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_16_1.png
[12]:
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=cf1['sweep_1'].latitude.values,
                                     central_longitude=cf1['sweep_1'].longitude.values,
                                    )
fig = pl.figure(figsize=(10,8))
ax = fig.add_subplot(111, projection=map_proj)
ax.set_boundary(circle, transform=ax.transAxes)

pm = cf1['sweep_1'].DBZH.wradlib.plot_ppi(proj=map_proj, ax=ax)
ax = pl.gca()
ax.gridlines(crs=map_proj)
[12]:
<cartopy.mpl.gridliner.Gridliner at 0x7f49a058ff40>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_17_1.png
[13]:
fig = pl.figure(figsize=(10, 8))
proj=ccrs.AzimuthalEquidistant(central_latitude=cf1['sweep_1'].latitude.values,
                               central_longitude=cf1['sweep_1'].longitude.values)
ax = fig.add_subplot(111, projection=proj)
pm = cf1['sweep_1'].DBZH.wradlib.plot_ppi(ax=ax)
ax.gridlines()
[13]:
<cartopy.mpl.gridliner.Gridliner at 0x7f49a3e34d60>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_18_1.png
[14]:
dbz = cf1['sweep_1']
dbz.DBZH.wradlib.plot_ppi()
[14]:
<matplotlib.collections.QuadMesh at 0x7f49a3fa3b20>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_19_1.png

Inspect radar moments

The dataarrays can be accessed by key or by attribute. Each dataarray has the datasets dimensions and coordinates of it’s parent dataset. There are attributes connected which are defined by Cf/Radial and/or ODIM_H5 standard.

[15]:
cf1['sweep_1'].DBZH
[15]:
<xarray.DataArray 'DBZH' (time: 360, range: 320)>
[115200 values with dtype=float32]
Coordinates:
    sweep_mode  <U20 ...
    latitude    float32 ...
    altitude    float32 ...
    longitude   float32 ...
    elevation   (time) float32 ...
    azimuth     (time) float32 ...
  * range       (range) float32 500.0 1500.0 2500.0 ... 318500.0 319500.0
    y           (time, range) float32 ...
    z           (time, range) float32 ...
    gr          (time, range) float32 ...
    rays        (time, range) float32 ...
    bins        (time, range) float32 ...
    x           (time, range) float32 ...
  * time        (time) datetime64[ns] 2011-06-10T11:40:06.694446592 ... 2011-06-10T11:40:06.638891008
Attributes:
    IMAGE_VERSION:  1.2
    standard_name:  radar_equivalent_reflectivity_factor_h
    long_name:      Equivalent reflectivity factor H
    units:          dBZ
[16]:
cf1['sweep_1']
[16]:
<xarray.Dataset>
Dimensions:       (range: 320, time: 360)
Coordinates:
    sweep_mode    <U20 ...
    latitude      float32 ...
    altitude      float32 ...
    longitude     float32 ...
    elevation     (time) float32 ...
    azimuth       (time) float32 ...
  * range         (range) float32 500.0 1500.0 2500.0 ... 318500.0 319500.0
    y             (time, range) float32 ...
    z             (time, range) float32 ...
    gr            (time, range) float32 ...
    rays          (time, range) float32 ...
    bins          (time, range) float32 ...
    x             (time, range) float32 ...
  * time          (time) datetime64[ns] 2011-06-10T11:40:06.694446592 ... 2011-06-10T11:40:06.638891008
Data variables:
    DBZH          (time, range) float32 ...
    sweep_number  int64 ...
    follow_mode   <U4 ...
    prt_mode      <U5 ...
    fixed_angle   float32 ...

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.

[17]:
cf1['sweep_1'].DBZH.sortby('time').plot(add_labels=False)
[17]:
<matplotlib.collections.QuadMesh at 0x7f49a3f555e0>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_24_1.png
[18]:
pm = cf1['sweep_1'].DBZH.wradlib.plot_ppi(proj={'latmin': 33e3})
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_25_0.png
[19]:
cf1.to_odim('knmi_odim.h5')
cf1.to_cfradial2('knmi_odim_as_cfradial.nc')

Import again

[20]:
cf1a = OdimH5('knmi_odim.h5', standard='cf', georef=True)
cf1b = CfRadial('knmi_odim_as_cfradial.nc', georef=True)
[21]:
cf1a['sweep_1']
[21]:
<xarray.Dataset>
Dimensions:       (range: 320, time: 360)
Coordinates:
    sweep_mode    <U20 ...
    latitude      float32 ...
    altitude      float32 ...
    longitude     float32 ...
    elevation     (time) float32 ...
    azimuth       (time) float32 ...
  * range         (range) float32 500.0 1500.0 2500.0 ... 318500.0 319500.0
    y             (time, range) float32 ...
    z             (time, range) float32 ...
    gr            (time, range) float32 ...
    rays          (time, range) float32 ...
    bins          (time, range) float32 ...
    x             (time, range) float32 ...
  * time          (time) datetime64[ns] 2011-06-10T11:40:06.694446592 ... 2011-06-10T11:40:06.638891008
Data variables:
    DBZH          (time, range) float32 ...
    sweep_number  int64 ...
    follow_mode   <U4 ...
    prt_mode      <U5 ...
    fixed_angle   float32 ...

Check equality

[22]:
xr.testing.assert_equal(cf1.root, cf1a.root)
xr.testing.assert_equal(cf1['sweep_1'], cf1a['sweep_1'])
xr.testing.assert_equal(cf1.root, cf1b.root)
xr.testing.assert_equal(cf1['sweep_1'], cf1b['sweep_1'])

Mask some values

[23]:
cf1['sweep_1']['DBZH'] = cf1['sweep_1']['DBZH'].where(cf1['sweep_1']['DBZH'] >= 0)
cf1['sweep_1']['DBZH'].sortby('time').plot()
[23]:
<matplotlib.collections.QuadMesh at 0x7f49a3cb71c0>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_33_1.png

Load Cf/Radial1 Volume Data

[24]:
fpath = 'netcdf/cfrad.20080604_002217_000_SPOL_v36_SUR.nc'
f = wrl.util.get_wradlib_data_file(fpath)
cf2 = CfRadial(f, georef=True)

Inspect root group

[25]:
cf2.root
[25]:
<xarray.Dataset>
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 ...
    platform_type        |S32 ...
    primary_axis         |S32 ...
    status_xml           |S1 ...
    instrument_type      |S32 ...
    time_coverage_start  |S32 ...
    time_coverage_end    |S32 ...
    latitude             float64 ...
    longitude            float64 ...
    altitude             float64 ...
    sweep_fixed_angle    (sweep) float32 ...
    sweep_group_name     (sweep) <U7 'sweep_1' 'sweep_2' ... 'sweep_8' 'sweep_9'
Attributes:
    Conventions:         CF/Radial instrument_parameters radar_parameters rad...
    version:             1.2
    title:               TIMREX
    institution:
    references:
    source:
    history:
    comment:
    instrument_name:     SPOLRVP8
    site_name:
    scan_name:
    scan_id:             0
    platform_is_mobile:  false
    n_gates_vary:        false

Inspect sweep group(s)

[26]:
cf2['sweep_1']
[26]:
<xarray.Dataset>
Dimensions:             (range: 996, time: 482)
Coordinates:
    sweep_mode          <U20 'azimuth_surveillance'
  * time                (time) datetime64[ns] 2008-06-04T00:15:03 ... 2008-06-04T00:15:50
  * range               (range) float32 150.0 300.0 ... 149250.0 149400.0
    azimuth             (time) float32 121.5 122.25 123.0 ... 121.5 122.25
    elevation           (time) float32 0.379 0.2362 0.1648 ... 0.5109 0.5109
    longitude           float64 120.4
    latitude            float64 22.53
    altitude            float64 45.0
    x                   (time, range) float32 127.89255 255.78506 ... 126313.25
    y                   (time, range) float32 -78.37265 -156.74529 ... -79697.73
    z                   (time, range) float32 45.0 46.0 47.0 ... 2685.0 2689.0
    gr                  (time, range) float32 149.99593 299.99182 ... 149354.5
    rays                (time, range) float32 121.5 121.5 ... 122.25 122.25
    bins                (time, range) float32 150.0 300.0 ... 149250.0 149400.0
Data variables:
    sweep_number        int32 ...
    polarization_mode   |S32 ...
    prt_mode            |S32 ...
    follow_mode         |S32 ...
    fixed_angle         float32 ...
    target_scan_rate    float32 ...
    pulse_width         (time) timedelta64[ns] ...
    prt                 (time) timedelta64[ns] ...
    nyquist_velocity    (time) float32 ...
    unambiguous_range   (time) float32 ...
    antenna_transition  (time) int8 ...
    n_samples           (time) int32 ...
    r_calib_index       (time) int8 ...
    scan_rate           (time) float32 ...
    DBZ                 (time, range) float32 ...
    VR                  (time, range) float32 ...

Inspect radar moments

[27]:
cf2['sweep_1'].DBZ
[27]:
<xarray.DataArray 'DBZ' (time: 482, range: 996)>
[480072 values with dtype=float32]
Coordinates:
    sweep_mode  <U20 'azimuth_surveillance'
  * time        (time) datetime64[ns] 2008-06-04T00:15:03 ... 2008-06-04T00:15:50
  * range       (range) float32 150.0 300.0 450.00003 ... 149250.0 149400.0
    azimuth     (time) float32 121.5 122.25 123.0 123.75 ... 120.75 121.5 122.25
    elevation   (time) float32 0.379 0.2362 0.1648 ... 0.5109 0.5109 0.5109
    longitude   float64 120.4
    latitude    float64 22.53
    altitude    float64 45.0
    x           (time, range) float32 127.89255 255.78506 ... 126313.25
    y           (time, range) float32 -78.37265 -156.74529 ... -79697.73
    z           (time, range) float32 45.0 46.0 47.0 ... 2681.0 2685.0 2689.0
    gr          (time, range) float32 149.99593 299.99182 ... 149204.6 149354.5
    rays        (time, range) float32 121.5 121.5 121.5 ... 122.25 122.25 122.25
    bins        (time, range) float32 150.0 300.0 ... 149250.0 149400.0
Attributes:
    long_name:             Computed Horizontal Co-polar Reflectivit
    standard_name:         equivalent_reflectivity_factor
    units:                 dBZ
    threshold_field_name:
    threshold_value:       -9999.0
    sampling_ratio:        1.0
    grid_mapping:          grid_mapping

Create simple plot

[28]:
cf2['sweep_1'].DBZ.plot()
[28]:
<matplotlib.collections.QuadMesh at 0x7f49a3aff520>
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_43_1.png
[29]:
cf2['sweep_1'].DBZ.plot.pcolormesh(x='x', y='y', add_labels=False)
pl.gca().set_aspect('equal')
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_44_0.png

Use wradlib DataArray connector

[30]:
pm = cf2['sweep_1'].DBZ.wradlib.plot_ppi()
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_46_0.png
[31]:
pm = cf2['sweep_1'].DBZ.wradlib.plot_ppi(proj='cg')
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_47_0.png

Export data to Cf/Radial2 and ODIM_H5

[32]:
cf2.to_cfradial2('timrex_cfradial2.nc')
cf2.to_odim('timrex_cfradial_as_odim.h5')

Import again

[33]:
cf2a = CfRadial('timrex_cfradial2.nc', georef=True)
cf2b = OdimH5('timrex_cfradial_as_odim.h5', standard='cf', georef=True)
[34]:
cf2a['sweep_1'].DBZ.plot.pcolormesh(x='x', y='y', add_labels=False)
pl.gca().set_aspect('equal')
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_52_0.png
[35]:
cf2b['sweep_1'].DBZ.plot.pcolormesh(x='x', y='y', add_labels=False)
pl.gca().set_aspect('equal')
../../_images/notebooks_fileio_wradlib_xarray_radial_odim_53_0.png

Check equality

For Cf/Radial there are issues with nan, which need to be fixed. For the ODIM_H5 intercomparison there are too problems with nan and issues with attributes.

[36]:
xr.testing.assert_equal(cf2.root, cf2a.root)
xr.testing.assert_equal(cf2['sweep_1'].drop(['DBZ', 'VR']),
                        cf2a['sweep_1'].drop(['DBZ', 'VR']))

xr.testing.assert_allclose(cf2.root.time_coverage_start,
                           cf2b.root.time_coverage_start)