CfRadial1#

[1]:
import xarray as xr
from open_radar_data import DATASETS

import xradar as xd

Download#

Fetching CfRadial1 radar data file from open-radar-data repository.

[2]:
filename = DATASETS.fetch("cfrad.20080604_002217_000_SPOL_v36_SUR.nc")
Downloading file 'cfrad.20080604_002217_000_SPOL_v36_SUR.nc' from 'https://github.com/openradar/open-radar-data/raw/main/data/cfrad.20080604_002217_000_SPOL_v36_SUR.nc' to '/home/docs/.cache/open-radar-data'.

xr.open_dataset#

Making use of the xarray cfradial1 backend. We also need to provide the group.

[3]:
ds = xr.open_dataset(filename, group="sweep_0", engine="cfradial1")
display(ds)
<xarray.Dataset> Size: 4MB
Dimensions:                    (azimuth: 483, range: 996)
Coordinates:
    time                       (azimuth) datetime64[ns] 4kB ...
  * range                      (range) float32 4kB 150.0 300.0 ... 1.494e+05
  * azimuth                    (azimuth) float32 2kB 0.0 0.75 ... 358.5 359.2
    elevation                  (azimuth) float32 2kB ...
    latitude                   float64 8B ...
    longitude                  float64 8B ...
    altitude                   float64 8B ...
Data variables: (12/18)
    sweep_number               int32 4B ...
    sweep_mode                 <U20 80B ...
    prt_mode                   |S32 32B ...
    follow_mode                |S32 32B ...
    sweep_fixed_angle          float32 4B ...
    pulse_width                (azimuth) timedelta64[ns] 4kB ...
    ...                         ...
    r_calib_index              (azimuth) int8 483B ...
    measured_transmit_power_h  (azimuth) float32 2kB ...
    measured_transmit_power_v  (azimuth) float32 2kB ...
    scan_rate                  (azimuth) float32 2kB ...
    DBZ                        (azimuth, range) float32 2MB ...
    VR                         (azimuth, range) float32 2MB ...

Plot Time vs. Azimuth#

We need to sort by time and specify the coordinate.

[4]:
ds.azimuth.plot(y="time")
[4]:
[<matplotlib.lines.Line2D at 0x7f9117cbecd0>]
../_images/notebooks_CfRadial1_7_1.png

Plot Range vs. Time#

We need to sort by time and specify the coordinate.

[5]:
ds.DBZ.sortby("time").plot(y="time")
[5]:
<matplotlib.collections.QuadMesh at 0x7f91167d2090>
../_images/notebooks_CfRadial1_9_1.png

Plot Range vs. Azimuth#

[6]:
ds.DBZ.plot()
[6]:
<matplotlib.collections.QuadMesh at 0x7f910d654cd0>
../_images/notebooks_CfRadial1_11_1.png

backend_kwargs#

The cfradial1 backend can be parameterized via kwargs. Please observe the possibilities below.

[7]:
?xd.io.CfRadial1BackendEntrypoint
[8]:
ds = xr.open_dataset(filename, group="sweep_0", engine="cfradial1", first_dim="time")
display(ds)
<xarray.Dataset> Size: 4MB
Dimensions:                    (time: 483, range: 996)
Coordinates:
  * time                       (time) datetime64[ns] 4kB 2008-06-04T00:15:03 ...
  * range                      (range) float32 4kB 150.0 300.0 ... 1.494e+05
    azimuth                    (time) float32 2kB ...
    elevation                  (time) float32 2kB ...
    latitude                   float64 8B ...
    longitude                  float64 8B ...
    altitude                   float64 8B ...
Data variables: (12/18)
    sweep_number               int32 4B ...
    sweep_mode                 <U20 80B ...
    prt_mode                   |S32 32B ...
    follow_mode                |S32 32B ...
    sweep_fixed_angle          float32 4B ...
    pulse_width                (time) timedelta64[ns] 4kB ...
    ...                         ...
    r_calib_index              (time) int8 483B ...
    measured_transmit_power_h  (time) float32 2kB ...
    measured_transmit_power_v  (time) float32 2kB ...
    scan_rate                  (time) float32 2kB ...
    DBZ                        (time, range) float32 2MB ...
    VR                         (time, range) float32 2MB ...
[9]:
ds = xr.open_dataset(
    filename, group="sweep_1", engine="cfradial1", backend_kwargs=dict(first_dim="time")
)
display(ds)
<xarray.Dataset> Size: 4MB
Dimensions:                    (time: 483, range: 996)
Coordinates:
  * time                       (time) datetime64[ns] 4kB 2008-06-04T00:15:50 ...
  * range                      (range) float32 4kB 150.0 300.0 ... 1.494e+05
    azimuth                    (time) float32 2kB ...
    elevation                  (time) float32 2kB ...
    latitude                   float64 8B ...
    longitude                  float64 8B ...
    altitude                   float64 8B ...
Data variables: (12/18)
    sweep_number               int32 4B ...
    sweep_mode                 <U20 80B ...
    prt_mode                   |S32 32B ...
    follow_mode                |S32 32B ...
    sweep_fixed_angle          float32 4B ...
    pulse_width                (time) timedelta64[ns] 4kB ...
    ...                         ...
    r_calib_index              (time) int8 483B ...
    measured_transmit_power_h  (time) float32 2kB ...
    measured_transmit_power_v  (time) float32 2kB ...
    scan_rate                  (time) float32 2kB ...
    DBZ                        (time, range) float32 2MB ...
    VR                         (time, range) float32 2MB ...

open_cfradial1_datatree#

The same works analoguous with the datatree loader. But additionally we can provide a sweep number or list.

[10]:
?xd.io.open_cfradial1_datatree
[11]:
dtree = xd.io.open_cfradial1_datatree(
    filename,
    first_dim="time",
    optional=False,
)
display(dtree)
<xarray.DatasetView> Size: 380B
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 4B ...
    time_coverage_start  |S32 32B ...
    time_coverage_end    |S32 32B ...
    latitude             float64 8B ...
    longitude            float64 8B ...
    altitude             float64 8B ...
    sweep_group_name     (sweep) <U7 252B 'sweep_0' 'sweep_1' ... 'sweep_8'
    sweep_fixed_angle    (sweep) float32 36B ...
Attributes:
    Conventions:         CF/Radial instrument_parameters radar_parameters rad...
    version:             1.2
    title:               TIMREX
    institution:
    references:
    source:
    history:
    comment:
    instrument_name:     SPOLRVP8
    platform_is_mobile:  false

Plot Sweep Range vs. Time#

[12]:
dtree["sweep_0"].ds.DBZ.plot()
[12]:
<matplotlib.collections.QuadMesh at 0x7f910d711210>
../_images/notebooks_CfRadial1_20_1.png

Plot Sweep Range vs. Azimuth#

[13]:
dtree["sweep_0"].ds.DBZ.sortby("azimuth").plot(y="azimuth")
[13]:
<matplotlib.collections.QuadMesh at 0x7f910c49e3d0>
../_images/notebooks_CfRadial1_22_1.png
[14]:
dtree = xd.io.open_cfradial1_datatree(filename, sweep=[0, 1, 8])
display(dtree)
<xarray.DatasetView> Size: 477B
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 4B ...
    platform_type        |S32 32B ...
    primary_axis         |S32 32B ...
    status_str           |S1 1B ...
    instrument_type      |S32 32B ...
    time_coverage_start  |S32 32B ...
    time_coverage_end    |S32 32B ...
    latitude             float64 8B ...
    longitude            float64 8B ...
    altitude             float64 8B ...
    sweep_group_name     (sweep) <U7 252B 'sweep_0' 'sweep_1' ... 'sweep_8'
    sweep_fixed_angle    (sweep) float32 36B ...
Attributes: (12/13)
    Conventions:         CF/Radial instrument_parameters radar_parameters rad...
    version:             1.2
    title:               TIMREX
    institution:
    references:
    source:
    ...                  ...
    comment:
    instrument_name:     SPOLRVP8
    site_name:
    scan_name:
    scan_id:             0
    platform_is_mobile:  false
[15]:
dtree = xd.io.open_cfradial1_datatree(filename, sweep=["sweep_0", "sweep_4", "sweep_8"])
display(dtree)
<xarray.DatasetView> Size: 477B
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 4B ...
    platform_type        |S32 32B ...
    primary_axis         |S32 32B ...
    status_str           |S1 1B ...
    instrument_type      |S32 32B ...
    time_coverage_start  |S32 32B ...
    time_coverage_end    |S32 32B ...
    latitude             float64 8B ...
    longitude            float64 8B ...
    altitude             float64 8B ...
    sweep_group_name     (sweep) <U7 252B 'sweep_0' 'sweep_1' ... 'sweep_8'
    sweep_fixed_angle    (sweep) float32 36B ...
Attributes: (12/13)
    Conventions:         CF/Radial instrument_parameters radar_parameters rad...
    version:             1.2
    title:               TIMREX
    institution:
    references:
    source:
    ...                  ...
    comment:
    instrument_name:     SPOLRVP8
    site_name:
    scan_name:
    scan_id:             0
    platform_is_mobile:  false