CfRadial1 to CfRadial2#

A data model transformation#

In this notebook we show how to transform the CfRadial1 Data model to a CfRadial2 representation.

We use some internal functions to show how xradar is working inside.

Within this notebook we reference to the CfRadial2.1 draft. As long as the FM301 WMO standard is not finalized we will rely on the drafts presented.

[1]:
import os

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

Open CfRadial1 file using xr.open_dataset#

Making use of the xarray netcdf4 backend. We get back all data and metadata in one single CfRadial1 Dataset.

[3]:
ds = xr.open_dataset(filename, engine="netcdf4")
with xr.set_options(
    display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000
):
    display(ds.load())
<xarray.Dataset> Size: 35MB
Dimensions:                           (sweep: 9, r_calib: 1, time: 4343,
                                       range: 996)
Coordinates:
  * time                              (time) datetime64[ns] 35kB 2008-06-04T0...
  * range                             (range) float32 4kB 150.0 ... 1.494e+05
Dimensions without coordinates: sweep, r_calib
Data variables:
    volume_number                     int32 4B 36
    platform_type                     |S32 32B b'fixed'
    primary_axis                      |S32 32B b'axis_z'
    status_xml                        |S1 1B b''
    instrument_type                   |S32 32B b'radar'
    radar_antenna_gain_h              float32 4B 45.15
    radar_antenna_gain_v              float32 4B 45.15
    radar_beam_width_h                float32 4B 0.92
    radar_beam_width_v                float32 4B 0.92
    radar_rx_bandwidth                float32 4B -9.999e+03
    time_coverage_start               |S32 32B b'2008-06-04T00:15:03Z'
    time_coverage_end                 |S32 32B b'2008-06-04T00:22:17Z'
    azimuth_correction                float32 4B 0.0
    elevation_correction              float32 4B 0.0
    range_correction                  float32 4B 0.0
    longitude_correction              float32 4B 0.0
    latitude_correction               float32 4B 0.0
    pressure_altitude_correction      float32 4B 0.0
    altitude_correction               float32 4B 0.0
    eastward_velocity_correction      float32 4B 0.0
    northward_velocity_correction     float32 4B 0.0
    vertical_velocity_correction      float32 4B 0.0
    heading_correction                float32 4B 0.0
    roll_correction                   float32 4B 0.0
    pitch_correction                  float32 4B 0.0
    drift_correction                  float32 4B 0.0
    rotation_correction               float32 4B 0.0
    tilt_correction                   float32 4B 0.0
    latitude                          float64 8B 22.53
    longitude                         float64 8B 120.4
    altitude                          float64 8B 45.0
    sweep_number                      (sweep) int32 36B 0 1 2 3 4 5 6 7 8
    sweep_mode                        (sweep) |S32 288B b'azimuth_surveillanc...
    polarization_mode                 (sweep) |S32 288B b'not_set' ... b'not_...
    prt_mode                          (sweep) |S32 288B b'not_set' ... b'not_...
    follow_mode                       (sweep) |S32 288B b'not_set' ... b'not_...
    fixed_angle                       (sweep) float32 36B 0.4999 1.099 ... 12.8
    target_scan_rate                  (sweep) float32 36B -9.999e+03 ... -9.9...
    sweep_start_ray_index             (sweep) int32 36B 0 483 966 ... 3376 3860
    sweep_end_ray_index               (sweep) int32 36B 482 965 ... 3859 4342
    r_calib_time                      (r_calib) |S32 32B b'2008-06-04T00:15:03Z'
    r_calib_pulse_width               (r_calib) timedelta64[ns] 8B 00:00:00
    r_calib_xmit_power_h              (r_calib) float32 4B 145.3
    r_calib_xmit_power_v              (r_calib) float32 4B 145.3
    r_calib_two_way_waveguide_loss_h  (r_calib) float32 4B -36.89
    r_calib_two_way_waveguide_loss_v  (r_calib) float32 4B -36.89
    r_calib_two_way_radome_loss_h     (r_calib) float32 4B -9.999e+03
    r_calib_two_way_radome_loss_v     (r_calib) float32 4B -9.999e+03
    r_calib_receiver_mismatch_loss    (r_calib) float32 4B -9.999e+03
    r_calib_radar_constant_h          (r_calib) float32 4B -70.72
    r_calib_radar_constant_v          (r_calib) float32 4B -70.72
    r_calib_antenna_gain_h            (r_calib) float32 4B -9.999e+03
    r_calib_antenna_gain_v            (r_calib) float32 4B -9.999e+03
    r_calib_noise_hc                  (r_calib) float32 4B -114.7
    r_calib_noise_vc                  (r_calib) float32 4B -114.7
    r_calib_noise_hx                  (r_calib) float32 4B -9.999e+03
    r_calib_noise_vx                  (r_calib) float32 4B -9.999e+03
    r_calib_receiver_gain_hc          (r_calib) float32 4B 36.89
    r_calib_receiver_gain_vc          (r_calib) float32 4B 36.89
    r_calib_receiver_gain_hx          (r_calib) float32 4B -9.999e+03
    r_calib_receiver_gain_vx          (r_calib) float32 4B -9.999e+03
    r_calib_base_dbz_1km_hc           (r_calib) float32 4B -9.999e+03
    r_calib_base_dbz_1km_vc           (r_calib) float32 4B -9.999e+03
    r_calib_base_dbz_1km_hx           (r_calib) float32 4B -9.999e+03
    r_calib_base_dbz_1km_vx           (r_calib) float32 4B -9.999e+03
    r_calib_sun_power_hc              (r_calib) float32 4B -9.999e+03
    r_calib_sun_power_vc              (r_calib) float32 4B -9.999e+03
    r_calib_sun_power_hx              (r_calib) float32 4B -9.999e+03
    r_calib_sun_power_vx              (r_calib) float32 4B -9.999e+03
    r_calib_noise_source_power_h      (r_calib) float32 4B -9.999e+03
    r_calib_noise_source_power_v      (r_calib) float32 4B -9.999e+03
    r_calib_power_measure_loss_h      (r_calib) float32 4B -9.999e+03
    r_calib_power_measure_loss_v      (r_calib) float32 4B -9.999e+03
    r_calib_coupler_forward_loss_h    (r_calib) float32 4B -9.999e+03
    r_calib_coupler_forward_loss_v    (r_calib) float32 4B -9.999e+03
    r_calib_zdr_correction            (r_calib) float32 4B -9.999e+03
    r_calib_ldr_correction_h          (r_calib) float32 4B -9.999e+03
    r_calib_ldr_correction_v          (r_calib) float32 4B -9.999e+03
    r_calib_system_phidp              (r_calib) float32 4B -9.999e+03
    r_calib_test_power_h              (r_calib) float32 4B -9.999e+03
    r_calib_test_power_v              (r_calib) float32 4B -9.999e+03
    r_calib_receiver_slope_hc         (r_calib) float32 4B -9.999e+03
    r_calib_receiver_slope_vc         (r_calib) float32 4B -9.999e+03
    r_calib_receiver_slope_hx         (r_calib) float32 4B -9.999e+03
    r_calib_receiver_slope_vx         (r_calib) float32 4B -9.999e+03
    azimuth                           (time) float32 17kB 121.5 122.2 ... 213.0
    elevation                         (time) float32 17kB 0.379 0.2362 ... 12.81
    pulse_width                       (time) timedelta64[ns] 35kB 00:00:00 .....
    prt                               (time) timedelta64[ns] 35kB -1 days +21...
    prt_ratio                         (time) timedelta64[ns] 35kB -1 days +21...
    nyquist_velocity                  (time) float32 17kB 26.92 26.92 ... 26.92
    unambiguous_range                 (time) float32 17kB 1.5e+05 ... 1.5e+05
    antenna_transition                (time) int8 4kB 0 0 0 0 0 0 ... 0 0 0 0 0
    n_samples                         (time) int32 17kB 192 192 192 ... 192 192
    r_calib_index                     (time) int8 4kB -1 -1 -1 -1 ... -1 -1 -1
    measured_transmit_power_h         (time) float32 17kB -9.999e+03 ... -9.9...
    measured_transmit_power_v         (time) float32 17kB -9.999e+03 ... -9.9...
    scan_rate                         (time) float32 17kB -3.277e+04 ... -3.2...
    DBZ                               (time, range) float32 17MB -3.07 ... nan
    VR                                (time, range) float32 17MB -25.24 ... nan
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

Extract CfRadial2 Groups and Subgroups#

Now as we have the CfRadial1 Dataset we can work towards extracting the CfRadial2 groups and subgroups.

Extract CfRadial2 Root-Group#

The following sections present the details of the information in the top-level (root) group of the data set.

We use a convenience function to extract the CfRadial2 root group from the CfRadial1 Dataset. We can call this function with one kwarg:

  • optional=False - only mandatory data and metadata is imported, defaults to True

optional=True#

[4]:
root = xd.io.backends.cfradial1._get_required_root_dataset(ds)
with xr.set_options(
    display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000
):
    display(root.load())
<xarray.Dataset> Size: 477B
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 4B 36
    platform_type        |S32 32B b'fixed'
    primary_axis         |S32 32B b'axis_z'
    status_str           |S1 1B b''
    instrument_type      |S32 32B b'radar'
    time_coverage_start  |S32 32B b'2008-06-04T00:15:03Z'
    time_coverage_end    |S32 32B b'2008-06-04T00:22:17Z'
    latitude             float64 8B 22.53
    longitude            float64 8B 120.4
    altitude             float64 8B 45.0
    sweep_group_name     (sweep) <U7 252B 'sweep_0' 'sweep_1' ... 'sweep_8'
    sweep_fixed_angle    (sweep) float32 36B 0.4999 1.099 1.802 ... 9.102 12.8
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

optional=False#

[5]:
root = xd.io.backends.cfradial1._get_required_root_dataset(ds, optional=False)
with xr.set_options(
    display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000
):
    display(root)
<xarray.Dataset> Size: 380B
Dimensions:              (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int32 4B 36
    time_coverage_start  |S32 32B b'2008-06-04T00:15:03Z'
    time_coverage_end    |S32 32B b'2008-06-04T00:22:17Z'
    latitude             float64 8B 22.53
    longitude            float64 8B 120.4
    altitude             float64 8B 45.0
    sweep_group_name     (sweep) <U7 252B 'sweep_0' 'sweep_1' ... 'sweep_8'
    sweep_fixed_angle    (sweep) float32 36B 0.4999 1.099 1.802 ... 9.102 12.8
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

Extract Root-Group metadata groups#

The Cfradial2 Data Model has a notion of root group metadata groups. Those groups provide additional metadata covering other aspects of the radar system.

The radar_parameters sub-group#

This group holds radar parameters specific to a radar instrument. It’s implemented as dictionary where the value can be used to override the name.

[6]:
display(xd.model.radar_parameters_subgroup)
{'radar_antenna_gain_h': None,
 'radar_antenna_gain_v': None,
 'radar_beam_width_h': None,
 'radar_beam_width_v': None,
 'half_power_beam_width_h': 'radar_beam_width_h',
 'half_power_beam_width_v': 'radar_beam_width_v',
 'radar_receiver_bandwidth': None,
 'radar_rx_bandwidth': 'radar_receiver_bandwidth'}

Again we use a convenience function to extract the group.

[7]:
radar_parameters = xd.io.backends.cfradial1._get_subgroup(
    ds, xd.model.radar_parameters_subgroup
)
display(radar_parameters.load())
<xarray.Dataset> Size: 20B
Dimensions:                   ()
Data variables:
    radar_receiver_bandwidth  float32 4B -9.999e+03
    radar_antenna_gain_v      float32 4B 45.15
    radar_beam_width_h        float32 4B 0.92
    radar_antenna_gain_h      float32 4B 45.15
    radar_beam_width_v        float32 4B 0.92

The radar_calibration sub-group#

For a radar, a different calibration is required for each pulse width. Therefore the calibration variables are arrays. If only one calibration is available it is squeezed by the reader.

[8]:
display(xd.model.radar_calibration_subgroup)
{'calib_index': None,
 'time': None,
 'pulse_width': None,
 'antenna_gain_h': None,
 'antenna_gain_v': None,
 'ant_gain_h': 'antenna_gain_h',
 'ant_gain_v': 'antenna_gain_v',
 'xmit_power_h': None,
 'xmit_power_v': None,
 'tx_power_h': 'xmit_power_h',
 'tx_power_v': 'xmit_power_v',
 'two_way_waveguide_loss_h': None,
 'two_way_waveguide_loss_v': None,
 'two_way_radome_loss_h': None,
 'two_way_radome_loss_v': None,
 'receiver_mismatch_loss': None,
 'receiver_mismatch_loss_h': None,
 'receiver_mismatch_loss_v': None,
 'rx_loss_h': 'receiver_mismatch_loss_h',
 'rx_loss_v': 'receiver_mismatch_loss_v',
 'radar_constant_h': None,
 'radar_constant_v': None,
 'probert_jones_correction': None,
 'dielectric_factor_used': None,
 'noise_hc': None,
 'noise_vc': None,
 'noise_hx': None,
 'noise_vx': None,
 'receiver_gain_hc': None,
 'receiver_gain_vc': None,
 'receiver_gain_hx': None,
 'receiver_gain_vx': None,
 'base_1km_hc': None,
 'base_dbz_1km_hc': 'base_1km_hc',
 'base_1km_vc': None,
 'base_dbz_1km_vc': 'base_1km_vc',
 'base_1km_hx': None,
 'base_dbz_1km_hx': 'base_1km_hx',
 'base_1km_vx': None,
 'base_dbz_1km_vx': 'base_1km_vx',
 'sun_power_hc': None,
 'sun_power_vc': None,
 'sun_power_hx': None,
 'sun_power_vx': None,
 'noise_source_power_h': None,
 'noise_source_power_v': None,
 'noise_power_short_pulse_h': 'noise_source_power_h',
 'noise_power_short_pulse_v': 'noise_source_power_v',
 'noise_power_h': 'noise_source_power_h',
 'noise_power_v': 'noise_source_power_v',
 'power_measure_loss_h': None,
 'power_measure_loss_v': None,
 'coupler_forward_loss_h': None,
 'coupler_forward_loss_v': None,
 'zdr_correction': None,
 'ldr_correction_h': None,
 'ldr_correction_v': None,
 'system_phidp': None,
 'test_power_h': None,
 'test_power_v': None,
 'receiver_slope_hc': None,
 'receiver_slope_vc': None,
 'receiver_slope_hx': None,
 'receiver_slope_vx': None,
 'k_squared_water': None,
 'i0_dbm_hc': None,
 'i0_dbm_vc': None,
 'i0_dbm_hx': None,
 'i0_dbm_vx': None,
 'dynamic_range_db_hc': None,
 'dynamic_range_db_vc': None,
 'dynamic_range_db_hx': None,
 'dynamic_range_db_vx': None,
 'dbz_correction': None}

Again we use a convenience function to extract the group.

[9]:
radar_calibration = xd.io.backends.cfradial1._get_radar_calibration(ds)
with xr.set_options(display_expand_data_vars=True):
    if radar_calibration:
        display(radar_calibration.load())
<xarray.Dataset> Size: 212B
Dimensions:                   ()
Data variables: (12/45)
    time                      |S32 32B b'2008-06-04T00:15:03Z'
    pulse_width               timedelta64[ns] 8B 00:00:00
    xmit_power_h              float32 4B 145.3
    xmit_power_v              float32 4B 145.3
    two_way_waveguide_loss_h  float32 4B -36.89
    two_way_waveguide_loss_v  float32 4B -36.89
    ...                        ...
    test_power_h              float32 4B -9.999e+03
    test_power_v              float32 4B -9.999e+03
    receiver_slope_hc         float32 4B -9.999e+03
    receiver_slope_vc         float32 4B -9.999e+03
    receiver_slope_hx         float32 4B -9.999e+03
    receiver_slope_vx         float32 4B -9.999e+03

The georeference_correction sub-group#

The following additional variables are used to quantify errors in the georeference data for moving platforms. These are constant for a volume.

[10]:
display(xd.model.georeferencing_correction_subgroup)
{'azimuth_correction': None,
 'elevation_correction': None,
 'range_correction': None,
 'longitude_correction': None,
 'latitude_correction': None,
 'pressure_altitude_correction': None,
 'altitude_correction': 'radar_altitude_correction',
 'radar_altitude_correction': None,
 'eastward_velocity_correction': 'eastward_ground_speed_correction',
 'eastward_ground_speed_correction': None,
 'northward_velocity_correction': 'northward_ground_speed_correction',
 'northward_ground_speed_correction': None,
 'vertical_velocity_correction': None,
 'heading_correction': None,
 'roll_correction': None,
 'pitch_correction': None,
 'drift_correction': None,
 'rotation_correction': None,
 'tilt_correction': None}

Again we use a convenience function to extract the group.

[11]:
georeference_correction = xd.io.backends.cfradial1._get_subgroup(
    ds, xd.model.georeferencing_correction_subgroup
)
with xr.set_options(display_expand_data_vars=True):
    display(georeference_correction.load())
<xarray.Dataset> Size: 64B
Dimensions:                            ()
Data variables: (12/16)
    heading_correction                 float32 4B 0.0
    drift_correction                   float32 4B 0.0
    longitude_correction               float32 4B 0.0
    vertical_velocity_correction       float32 4B 0.0
    tilt_correction                    float32 4B 0.0
    northward_ground_speed_correction  float32 4B 0.0
    ...                                 ...
    radar_altitude_correction          float32 4B 0.0
    latitude_correction                float32 4B 0.0
    roll_correction                    float32 4B 0.0
    azimuth_correction                 float32 4B 0.0
    range_correction                   float32 4B 0.0
    pressure_altitude_correction       float32 4B 0.0

Sweep groups#

This section provides details of the information in each sweep group. The name of the sweep groups is found in the sweep_group_name array variable in the root group.

[12]:
root.sweep_group_name
[12]:
<xarray.DataArray 'sweep_group_name' (sweep: 9)> Size: 252B
array(['sweep_0', 'sweep_1', 'sweep_2', 'sweep_3', 'sweep_4', 'sweep_5',
       'sweep_6', 'sweep_7', 'sweep_8'], dtype='<U7')
Dimensions without coordinates: sweep

Again we use a convenience function to extract the different sweep groups. We can call this function with kwargs:

  • optional=False - only mandatory data and metadata is imported, defaults to True

  • first_dim="time - return first dimension as time, defaults toauto (return either as azimuth (PPI) or elevation (RHI)to time

  • site_coords=False - do not add radar site coordinates to the Sweep-Dataset, defaults to True

Examining first sweep with default kwargs.#

[13]:
sweeps = xd.io.backends.cfradial1._get_sweep_groups(ds)
with xr.set_options(display_expand_data_vars=True):
    display(sweeps["sweep_0"])
<xarray.Dataset> Size: 4MB
Dimensions:                    (azimuth: 483, range: 996)
Coordinates:
    time                       (azimuth) datetime64[ns] 4kB 2008-06-04T00:15:...
  * 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 0.5164 0.5219 ... 0.5219
    latitude                   float64 8B 22.53
    longitude                  float64 8B 120.4
    altitude                   float64 8B 45.0
Data variables: (12/18)
    sweep_number               int32 4B 0
    sweep_mode                 <U20 80B 'azimuth_surveillance'
    prt_mode                   |S32 32B b'not_set'
    follow_mode                |S32 32B b'not_set'
    sweep_fixed_angle          float32 4B 0.4999
    pulse_width                (azimuth) timedelta64[ns] 4kB 00:00:00 ... 00:...
    ...                         ...
    r_calib_index              (azimuth) int8 483B -1 -1 -1 -1 ... -1 -1 -1 -1
    measured_transmit_power_h  (azimuth) float32 2kB -9.999e+03 ... -9.999e+03
    measured_transmit_power_v  (azimuth) float32 2kB -9.999e+03 ... -9.999e+03
    scan_rate                  (azimuth) float32 2kB -3.277e+04 ... -3.277e+04
    DBZ                        (azimuth, range) float32 2MB 20.7 39.97 ... -1.32
    VR                         (azimuth, range) float32 2MB -23.38 ... -25.39

Examining first sweep with optional=False#

[14]:
sweeps = xd.io.backends.cfradial1._get_sweep_groups(ds, optional=False)
with xr.set_options(display_expand_data_vars=True):
    display(sweeps["sweep_0"])
<xarray.Dataset> Size: 4MB
Dimensions:            (azimuth: 483, range: 996)
Coordinates:
    time               (azimuth) datetime64[ns] 4kB 2008-06-04T00:15:34 ... 2...
  * range              (range) float32 4kB 150.0 300.0 ... 1.492e+05 1.494e+05
  * azimuth            (azimuth) float32 2kB 0.0 0.75 1.5 ... 357.8 358.5 359.2
    elevation          (azimuth) float32 2kB 0.5164 0.5219 ... 0.5219 0.5219
    latitude           float64 8B 22.53
    longitude          float64 8B 120.4
    altitude           float64 8B 45.0
Data variables:
    sweep_number       int32 4B 0
    sweep_mode         <U20 80B 'azimuth_surveillance'
    prt_mode           |S32 32B b'not_set'
    follow_mode        |S32 32B b'not_set'
    sweep_fixed_angle  float32 4B 0.4999
    DBZ                (azimuth, range) float32 2MB 20.7 39.97 ... -3.16 -1.32
    VR                 (azimuth, range) float32 2MB -23.38 10.22 ... -25.39

optional=False and site_coords=False#

[15]:
sweeps = xd.io.backends.cfradial1._get_sweep_groups(
    ds, optional=False, site_coords=False
)
with xr.set_options(display_expand_data_vars=True):
    display(sweeps["sweep_0"])
<xarray.Dataset> Size: 4MB
Dimensions:            (azimuth: 483, range: 996)
Coordinates:
    time               (azimuth) datetime64[ns] 4kB 2008-06-04T00:15:34 ... 2...
  * range              (range) float32 4kB 150.0 300.0 ... 1.492e+05 1.494e+05
  * azimuth            (azimuth) float32 2kB 0.0 0.75 1.5 ... 357.8 358.5 359.2
    elevation          (azimuth) float32 2kB 0.5164 0.5219 ... 0.5219 0.5219
Data variables:
    sweep_number       int32 4B 0
    sweep_mode         <U20 80B 'azimuth_surveillance'
    prt_mode           |S32 32B b'not_set'
    follow_mode        |S32 32B b'not_set'
    sweep_fixed_angle  float32 4B 0.4999
    DBZ                (azimuth, range) float32 2MB 20.7 39.97 ... -3.16 -1.32
    VR                 (azimuth, range) float32 2MB -23.38 10.22 ... -25.39

optional=False, site_coords=True and first_dim="auto"#

[16]:
sweeps = xd.io.backends.cfradial1._get_sweep_groups(
    ds, optional=False, site_coords=False, first_dim="time"
)
with xr.set_options(display_expand_data_vars=True):
    display(sweeps["sweep_0"])
<xarray.Dataset> Size: 4MB
Dimensions:            (time: 483, range: 996)
Coordinates:
  * time               (time) datetime64[ns] 4kB 2008-06-04T00:15:03 ... 2008...
  * range              (range) float32 4kB 150.0 300.0 ... 1.492e+05 1.494e+05
    azimuth            (time) float32 2kB 121.5 122.2 123.0 ... 122.2 123.0
    elevation          (time) float32 2kB 0.379 0.2362 0.1648 ... 0.5109 0.5109
Data variables:
    sweep_number       int32 4B 0
    sweep_mode         <U20 80B 'azimuth_surveillance'
    prt_mode           |S32 32B b'not_set'
    follow_mode        |S32 32B b'not_set'
    sweep_fixed_angle  float32 4B 0.4999
    DBZ                (time, range) float32 2MB -3.07 16.3 ... -0.2697 -0.01028
    VR                 (time, range) float32 2MB -25.24 3.48 ... 3.19 -23.05

Read as CfRadial2 data representation#

xradar provides two easy ways to retrieve the CfRadial1 data as CfRadial2 groups.

DataTree#

This is the most complete representation as a DataTree. All groups and subgroups are represented in a tree-like structure. Can be parameterized using kwargs. Easy write to netCDF4.

[17]:
dtree = xd.io.open_cfradial1_datatree(filename)
with xr.set_options(display_expand_data_vars=True, display_expand_attrs=True):
    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

Each DataTree-node itself represents another DataTree.

[18]:
display(dtree["radar_parameters"].load())
<xarray.DatasetView> Size: 20B
Dimensions:                   (sweep: 9)
Dimensions without coordinates: sweep
Data variables:
    radar_receiver_bandwidth  float32 4B -9.999e+03
    radar_antenna_gain_v      float32 4B 45.15
    radar_beam_width_h        float32 4B 0.92
    radar_antenna_gain_h      float32 4B 45.15
    radar_beam_width_v        float32 4B 0.92
[19]:
with xr.set_options(display_expand_data_vars=True):
    display(dtree["sweep_7"].load())
<xarray.DatasetView> Size: 4MB
Dimensions:                    (sweep: 9, azimuth: 484, range: 996)
Coordinates:
    time                       (azimuth) datetime64[ns] 4kB 2008-06-04T00:21:...
  * 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 9.102 9.102 ... 9.102 9.102
    latitude                   float64 8B 22.53
    longitude                  float64 8B 120.4
    altitude                   float64 8B 45.0
Dimensions without coordinates: sweep
Data variables: (12/18)
    sweep_number               int32 4B 7
    sweep_mode                 <U20 80B 'azimuth_surveillance'
    prt_mode                   |S32 32B b'not_set'
    follow_mode                |S32 32B b'not_set'
    sweep_fixed_angle          float32 4B 9.102
    pulse_width                (azimuth) timedelta64[ns] 4kB 00:00:00 ... 00:...
    ...                         ...
    r_calib_index              (azimuth) int8 484B -1 -1 -1 -1 ... -1 -1 -1 -1
    measured_transmit_power_h  (azimuth) float32 2kB -9.999e+03 ... -9.999e+03
    measured_transmit_power_v  (azimuth) float32 2kB -9.999e+03 ... -9.999e+03
    scan_rate                  (azimuth) float32 2kB -3.277e+04 ... -3.277e+04
    DBZ                        (azimuth, range) float32 2MB 21.88 ... -8.661
    VR                         (azimuth, range) float32 2MB -0.4204 ... 11.4