Supported radar data formats

The binary encoding of many radar products is a major obstacle for many potential radar users. Often, decoder software is not easily available. In case formats are documented, the implementation of decoders is a major programming effort. This tutorial provides an overview of the data formats currently supported by \(\omega radlib\). We seek to continuously enhance the range of supported formats, so this document is only a snapshot. If you need a specific file format to be supported by \(\omega radlib\), please raise an issue of type enhancement. You can provide support by adding documents which help to decode the format, e.g. format reference documents or software code in other languages for decoding the format.

At the moment, supported format means that the radar format can be read and further processed by wradlib. Normally, wradlib will return a numpy array of data values and a dictionary of metadata - if the file contains any.

Warning

Since \(\omega radlib\) version 1.18 the xarray backend engines for polar radar data have been renamed and prepended with wradlib- (eg. cfradial1 -> wradlib-cfradial1). 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.

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.

In the following, we will provide an overview of file formats which can be currently read by \(\omega radlib\).

Reading weather radar files is done via the wradlib.io module. There you will find a complete function reference.

German Weather Service: DX format

The German Weather Service uses the DX file format to encode local radar sweeps. DX data are in polar coordinates. The naming convention is as follows:

raa00-dx_<location-id>-<YYMMDDHHMM>-<location-abreviation>---bin

or

raa00-dx_<location-id>-<YYYYMMDDHHMM>-<location-abreviation>---bin

Read and plot DX radar data from DWD provides an extensive introduction into working with DX data. For now, we would just like to know how to read the data:

[2]:
fpath = "dx/raa00-dx_10908-0806021655-fbg---bin.gz"
f = wrl.util.get_wradlib_data_file(fpath)
data, metadata = wrl.io.read_dx(f)

Here, data is a two dimensional array of shape (number of azimuth angles, number of range gates). This means that the number of rows of the array corresponds to the number of azimuth angles of the radar sweep while the number of columns corresponds to the number of range gates per ray.

[3]:
print(data.shape)
print(metadata.keys())
(360, 128)
dict_keys(['producttype', 'datetime', 'radarid', 'bytes', 'version', 'cluttermap', 'dopplerfilter', 'statfilter', 'elevprofile', 'message', 'elev', 'azim', 'clutter'])
[4]:
fig = pl.figure(figsize=(10, 10))
ax, im = wrl.vis.plot_ppi(data, fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_10_0.png

German Weather Service: RADOLAN (quantitative) composit

The quantitative composite format of the DWD (German Weather Service) was established in the course of the RADOLAN project. Most quantitative composite products from the DWD are distributed in this format, e.g. the R-series (RX, RY, RH, RW, …), the S-series (SQ, SH, SF, …), and the E-series (European quantitative composite, e.g. EZ, EH, EB). Please see the composite format description for a full reference and a full table of products (unfortunately only in German language). An extensive section covering many RADOLAN aspects is here: RADOLAN

Currently, the RADOLAN composites have a spatial resolution of 1km x 1km, with the national composits (R- and S-series) being 900 x 900 grids, and the European composits 1500 x 1400 grids. The projection is polar-stereographic. The products can be read by the following function:

Note

Since \(\omega radlib\) version 1.10 a RADOLAN reader wradlib.io.radolan.open_radolan_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_radolan_backend.

[5]:
fpath = "radolan/misc/raa01-rw_10000-1408102050-dwd---bin.gz"
f = wrl.util.get_wradlib_data_file(fpath)
data, metadata = wrl.io.read_radolan_composite(f)
Downloading file 'radolan/misc/raa01-rw_10000-1408102050-dwd---bin.gz' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/radolan/misc/raa01-rw_10000-1408102050-dwd---bin.gz' to '/home/runner/work/wradlib-notebooks/wradlib-notebooks/wradlib-data'.

Here, data is a two dimensional integer array of shape (number of rows, number of columns). Different product types might need different levels of postprocessing, e.g. if the product contains rain rates or accumulations, you will normally have to divide data by factor 10. metadata is again a dictionary which provides metadata from the files header section, e.g. using the keys producttype, datetime, intervalseconds, nodataflag.

[6]:
print(data.shape)
print(metadata.keys())
(900, 900)
dict_keys(['producttype', 'datetime', 'radarid', 'datasize', 'formatversion', 'maxrange', 'radolanversion', 'precision', 'intervalseconds', 'nrow', 'ncol', 'radarlocations', 'nodataflag', 'secondary', 'nodatamask', 'cluttermask'])

Masking the NoData (or missing) values can be done by:

[7]:
maskeddata = np.ma.masked_equal(data, metadata["nodataflag"])
[8]:
fig = pl.figure(figsize=(10, 8))
# get coordinates
radolan_grid_xy = wrl.georef.get_radolan_grid(900, 900)
x = radolan_grid_xy[:, :, 0]
y = radolan_grid_xy[:, :, 1]

# create quick plot with colorbar and title
pl.figure(figsize=(10, 8))
pl.pcolormesh(x, y, maskeddata)
[8]:
<matplotlib.collections.QuadMesh at 0x7f11bdd90910>
<Figure size 1000x800 with 0 Axes>
../../_images/notebooks_fileio_wradlib_radar_formats_18_2.png

HDF5

OPERA HDF5 (ODIM_H5)

HDF5 is a data model, library, and file format for storing and managing data. The OPERA program developed a convention (or information model) on how to store and exchange radar data in hdf5 format. It is based on the work of COST Action 717 and is used e.g. in real-time operations in the Nordic European countries. The OPERA Data and Information Model (ODIM) is documented e.g. in this report. Make use of these documents in order to understand the organization of OPERA hdf5 files!

Note

Since \(\omega radlib\) version 1.10 an ODIM_H5 reader wradlib.io.open_odim_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_odim_backend.

From \(\omega radlib\) version 1.19 ODIM_H5 reading code is imported from xradar-package whenever and wherever necessary.

Former xarray-based implementations will be deprecated in future versions.

The hierarchical nature of HDF5 can be described as being similar to directories, files, and links on a hard-drive. Actual metadata are stored as so-called attributes, and these attributes are organized together in so-called groups. Binary data are stored as so-called datasets. As for ODIM_H5, the root (or top level) group contains three groups of metadata: these are called what (object, information model version, and date/time information), where (geographical information), and how (quality and optional/recommended metadata). For a very simple product, e.g. a CAPPI, the data is organized in a group called dataset1 which contains another group called data1 where the actual binary data are found in data. In analogy with a file system on a hard-disk, the HDF5 file containing this simple product is organized like this:

/
/what
/where
/how
/dataset1
/dataset1/data1
/dataset1/data1/data

The philosophy behind the \(\omega radlib\) interface to OPERA’s data model is very straightforward: \(\omega radlib\) simply translates the complete file structure to one dictionary and returns this dictionary to the user. Thus, the potential complexity of the stored data is kept and it is left to the user how to proceed with this data. The keys of the output dictionary are strings that correspond to the “directory trees” shown above. Each key ending with /data points to a Dataset (i.e. a numpy array of data). Each key ending with /what, /where or /how points to another dictionary of metadata. The entire output can be obtained by:

[9]:
fpath = "hdf5/knmi_polar_volume.h5"
f = wrl.util.get_wradlib_data_file(fpath)
fcontent = wrl.io.read_opera_hdf5(f)

The user should inspect the output obtained from his or her hdf5 file in order to see how access those items which should be further processed. In order to get a readable overview of the output dictionary, one can use the pretty printing module:

[10]:
# which keyswords can be used to access the content?
print(fcontent.keys())
# print the entire content including values of data and metadata
# (numpy arrays will not be entirely printed)
print(fcontent["dataset1/data1/data"])
dict_keys(['dataset1/data1/data', 'dataset1/data1/what', 'dataset1/what', 'dataset1/where', 'dataset10/data1/data', 'dataset10/data1/what', 'dataset10/what', 'dataset10/where', 'dataset11/data1/data', 'dataset11/data1/what', 'dataset11/what', 'dataset11/where', 'dataset12/data1/data', 'dataset12/data1/what', 'dataset12/what', 'dataset12/where', 'dataset13/data1/data', 'dataset13/data1/what', 'dataset13/what', 'dataset13/where', 'dataset14/data1/data', 'dataset14/data1/what', 'dataset14/what', 'dataset14/where', 'dataset2/data1/data', 'dataset2/data1/what', 'dataset2/what', 'dataset2/where', 'dataset3/data1/data', 'dataset3/data1/what', 'dataset3/what', 'dataset3/where', 'dataset4/data1/data', 'dataset4/data1/what', 'dataset4/what', 'dataset4/where', 'dataset5/data1/data', 'dataset5/data1/what', 'dataset5/what', 'dataset5/where', 'dataset6/data1/data', 'dataset6/data1/what', 'dataset6/what', 'dataset6/where', 'dataset7/data1/data', 'dataset7/data1/what', 'dataset7/what', 'dataset7/where', 'dataset8/data1/data', 'dataset8/data1/what', 'dataset8/what', 'dataset8/where', 'dataset9/data1/data', 'dataset9/data1/what', 'dataset9/what', 'dataset9/where', 'what', 'where'])
[[107  97  47 ...   0   0   0]
 [111 112  45 ...   0   0   0]
 [134 147  87 ...   0   0   0]
 ...
 [109  91  37 ...   0   0   0]
 [109  91  45 ...   0   0   0]
 [107 100  40 ...   0   0   0]]

Please note that in order to experiment with such datasets, you can download hdf5 sample data from the OPERA or use the example data provided with the wradlib-data repository.

[11]:
fig = pl.figure(figsize=(10, 10))
im = wrl.vis.plot_ppi(fcontent["dataset1/data1/data"], fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_26_0.png

GAMIC HDF5

GAMIC refers to the commercial GAMIC Enigma MURAN software which exports data in hdf5 format. The concept is quite similar to the above OPERA HDF5 (ODIM_H5) format.

Note

Since \(\omega radlib\) version 1.10 an GAMIC reader wradlib.io.hdf.open_gamic_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_gamic_backend.

From \(\omega radlib\) version 1.19 GAMIC reading code is imported from xradar-package whenever and wherever necessary.

Former xarray-based implementations will be deprecated in future versions.

Such a file (typical ending: .mvol) can be read by:

[12]:
fpath = "hdf5/2014-08-10--182000.ppi.mvol"
f = wrl.util.get_wradlib_data_file(fpath)
data, metadata = wrl.io.read_gamic_hdf5(f)
Downloading file 'hdf5/2014-08-10--182000.ppi.mvol' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/hdf5/2014-08-10--182000.ppi.mvol' to '/home/runner/work/wradlib-notebooks/wradlib-notebooks/wradlib-data'.

While metadata represents the usual dictionary of metadata, the data variable is a dictionary which might contain several numpy arrays with the keywords of the dictionary indicating different moments.

[13]:
print(metadata.keys())
print(metadata["VOL"])
print(metadata["SCAN0"].keys())
dict_keys(['SCAN0', 'VOL'])
{'Latitude': 50.73052, 'Longitude': 7.071663, 'Height': 99.5}
dict_keys(['PRF', 'angle_step', 'angle_sync', 'azi_start', 'azi_stop', 'bin_count', 'elevation', 'filter', 'half_resolution', 'output64', 'pulse_width', 'radar_wave_length', 'range', 'range_samples', 'range_start', 'range_step', 'ray_count', 'scan_speed', 'time_samples', 'unfolding', 'bin_range', 'zero_index', 'az', 'el', 'r', 'Time', 'max_range'])
[14]:
print(data["SCAN0"].keys())
print(data["SCAN0"]["PHIDP"].keys())
print(data["SCAN0"]["PHIDP"]["data"].shape)
dict_keys(['KDP', 'PHIDP', 'ZH', 'ZV', 'RHOHV', 'UH', 'UV', 'VH', 'VV', 'WH', 'WV', 'ZDR'])
dict_keys(['data', 'dyn_range_max', 'dyn_range_min'])
(360, 1000)
[15]:
fig = pl.figure(figsize=(10, 10))
im = wrl.vis.plot_ppi(data["SCAN0"]["ZH"]["data"], fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_33_0.png

Generic HDF5

This is a generic hdf5 reader, which will read any hdf5 structure.

[16]:
fpath = "hdf5/2014-08-10--182000.ppi.mvol"
f = wrl.util.get_wradlib_data_file(fpath)
fcontent = wrl.io.read_generic_hdf5(f)
[17]:
print(fcontent.keys())
dict_keys(['how', 'scan0/how', 'scan0/how/extended', 'scan0/moment_0', 'scan0/moment_1', 'scan0/moment_10', 'scan0/moment_11', 'scan0/moment_2', 'scan0/moment_3', 'scan0/moment_4', 'scan0/moment_5', 'scan0/moment_6', 'scan0/moment_7', 'scan0/moment_8', 'scan0/moment_9', 'scan0/ray_header', 'scan0/what', 'what', 'where'])
[18]:
print(fcontent["where"])
print(fcontent["how"])
print(fcontent["scan0/moment_3"].keys())
print(fcontent["scan0/moment_3"]["attrs"])
print(fcontent["scan0/moment_3"]["data"].shape)
{'attrs': {'height': 99.5, 'lat': 50.73052, 'lon': 7.071663}}
{'attrs': {'azimuth_beam': 1.0, 'elevation_beam': 1.0, 'host_name': 'radar.meteo.uni-bonn.de', 'sdp_name': 'ENIGMA III DUALPOL', 'site_name': '12345', 'software': 'MURAN', 'template_name': 'ppi_1p5deg'}}
dict_keys(['attrs', 'data'])
{'dyn_range_max': 95.5, 'dyn_range_min': -32.0, 'format': 'UV8', 'moment': 'UH'}
(360, 1000)
[19]:
fig = pl.figure(figsize=(10, 10))
im = wrl.vis.plot_ppi(fcontent["scan0/moment_3"]["data"], fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_39_0.png

NetCDF

The NetCDF format also claims to be self-describing. However, as for all such formats, the developers of netCDF also admit that “[…] the mere use of netCDF is not sufficient to make data self-describing and meaningful to both humans and machines […]” (see here. Different radar operators or data distributors will use different naming conventions and data hierarchies (i.e. “data models”) that the reading program might need to know about.

\(\omega radlib\) provides two solutions to address this challenge. The first one ignores the concept of data models and just pulls all data and metadata from a NetCDF file (wradlib.io.read_generic_netcdf(). The second is designed for a specific data model used by the EDGE software (wradlib.io.read_edge_netcdf()).

Note

Since \(\omega radlib\) version 1.10 a Cf/Radial1 reader wradlib.io.xarray.open_cfradial1_dataset() and a Cf/Radial2 reader wradlib.io.xarray.open_cfradial2_dataset() for CF versions 1.X and 2 based on Xarray are available. Please read the more indepth notebooks wradlib_cfradial1_backend and wradlib_cfradial2_backend.

From \(\omega radlib\) version 1.19 Cf/Radial reading code is imported from xradar-package whenever and wherever necessary.

Generic NetCDF reader (includes CfRadial)

\(\omega radlib\) provides a function that will virtually read any NetCDF file irrespective of the data model: wradlib.io.read_generic_netcdf(). It is built upon Python’s netcdf4 library. wradlib.io.read_generic_netcdf() will return only one object, a dictionary, that contains all the contents of the NetCDF file corresponding to the original file structure. This includes all the metadata, as well as the so called “dimensions” (describing the dimensions of the actual data arrays) and the “variables” which will contains the actual data. Users can use this dictionary at will in order to query data and metadata; however, they should make sure to consider the documentation of the corresponding data model. wradlib.io.read_generic_netcdf() has been shown to work with a lot of different data models, most notably CfRadial (see here for details). A typical call to wradlib.io.read_generic_netcdf() would look like:

[20]:
fpath = "netcdf/example_cfradial_ppi.nc"
f = wrl.util.get_wradlib_data_file(fpath)
outdict = wrl.io.read_generic_netcdf(f)
for key in outdict.keys():
    print(key)
Downloading file 'netcdf/example_cfradial_ppi.nc' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/netcdf/example_cfradial_ppi.nc' to '/home/runner/work/wradlib-notebooks/wradlib-notebooks/wradlib-data'.
comment
title
Conventions
source
version
references
instrument_name
institution
field_names
history
dimensions
variables

Please see this example notebook to get started.

EDGE NetCDF

EDGE is a commercial software for radar control and data analysis provided by the Enterprise Electronics Corporation. It allows for netCDF data export. The resulting files can be read by wradlib.io.read_generic_netcdf(), but \(\omega radlib\) also provides a specific function, wradlib.io.read_edge_netcdf() to return metadata and data as seperate objects:

[21]:
fpath = "netcdf/edge_netcdf.nc"
f = wrl.util.get_wradlib_data_file(fpath)
data, metadata = wrl.io.read_edge_netcdf(f)
print(data.shape)
print(metadata.keys())
Downloading file 'netcdf/edge_netcdf.nc' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/netcdf/edge_netcdf.nc' to '/home/runner/work/wradlib-notebooks/wradlib-notebooks/wradlib-data'.
(360, 240)
dict_keys(['TypeName', 'DataType', 'Latitude', 'Longitude', 'Height', 'FractionalTime', 'attributes', 'NyquistVelocity-unit', 'NyquistVelocity-value', 'vcp-unit', 'vcp-value', 'radarName-unit', 'radarName-value', 'ColorMap-unit', 'ColorMap-value', 'Elevation', 'ElevationUnits', 'MissingData', 'RangeFolded', 'RadarParameters', 'PRF-unit', 'PRF-value', 'PulseWidth-unit', 'PulseWidth-value', 'MaximumRange-unit', 'MaximumRange-value', 'ConversionPlugin', 'az', 'r', 'sitecoords', 'time', 'max_range'])

Gematronik Rainbow

Rainbow refers to the commercial RAINBOW®5 APPLICATION SOFTWARE which exports data in an XML flavour, which due to binary data blobs violates XML standard. Gematronik provided python code for implementing this reader in \(\omega radlib\), which is very much appreciated.

The philosophy behind the \(\omega radlib\) interface to Gematroniks data model is very straightforward: \(\omega radlib\) simply translates the complete xml file structure to one dictionary and returns this dictionary to the user. Thus, the potential complexity of the stored data is kept and it is left to the user how to proceed with this data. The keys of the output dictionary are strings that correspond to the “xml nodes” and “xml attributes”. Each data key points to a Dataset (i.e. a numpy array of data).

Note

Since \(\omega radlib\) version 1.13 a Rainbow5 reader wradlib.io.xarray.open_rainbow_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_rainbow_backend.

From \(\omega radlib\) version 1.19 Rainbow5 reading code is imported from xradar-package whenever and wherever necessary.

Such a file (typical ending: .vol or .azi) can be read by:

[22]:
fpath = "rainbow/2013070308340000dBuZ.azi"
f = wrl.util.get_wradlib_data_file(fpath)
fcontent = wrl.io.read_rainbow(f)

The user should inspect the output obtained from his or her Rainbow file in order to see how access those items which should be further processed. In order to get a readable overview of the output dictionary, one can use the pretty printing module:

[23]:
# which keyswords can be used to access the content?
print(fcontent.keys())
# print the entire content including values of data and metadata
# (numpy arrays will not be entirely printed)
print(fcontent["volume"]["sensorinfo"])
dict_keys(['volume'])
{'@type': 'rainscanner', '@id': 'WUE', '@name': 'Wuestebach', 'lon': '6.330970', 'lat': '50.504900', 'alt': '0.000000', 'wavelen': '0.05', 'beamwidth': '1'}

You can check this example notebook for getting a first impression.

Vaisala Sigmet IRIS

IRIS refers to the commercial Vaisala Sigmet Interactive Radar Information System. The Vaisala Sigmet Digital Receivers export data in a well documented binary format.

The philosophy behind the \(\omega radlib\) interface to the IRIS data model is very straightforward: \(\omega radlib\) simply translates the complete binary file structure to one dictionary and returns this dictionary to the user. Thus, the potential complexity of the stored data is kept and it is left to the user how to proceed with this data. The keys of the output dictionary are strings that correspond to the Sigmet Data Structures.

Note

Since \(\omega radlib\) version 1.12 an IRIS reader wradlib.io.iris.open_iris_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_iris_backend.

At the same time the reader has changed with respect to performance. So the ray metadata is only read once per sweep and is only included once in the output of read_iris. Currently we keep backwards compatibility, but this behaviour is deprecated and will be changed in a future version. See the two examples below.

From \(\omega radlib\) version 1.19 IRIS reading code is imported from xradar-package whenever and wherever necessary.

Such a file (typical ending: *.RAWXXXX) can be read by:

[24]:
fpath = "sigmet/cor-main131125105503.RAW2049"
f = wrl.util.get_wradlib_data_file(fpath)
fcontent = wrl.io.read_iris(f, keep_old_sweep_data=False)
Downloading file 'sigmet/cor-main131125105503.RAW2049' from 'https://github.com/wradlib/wradlib-data/raw/pooch/data/sigmet/cor-main131125105503.RAW2049' to '/home/runner/work/wradlib-notebooks/wradlib-notebooks/wradlib-data'.
/home/runner/micromamba-root/envs/wradlib-notebooks/lib/python3.11/site-packages/xradar/io/backends/iris.py:234: RuntimeWarning: invalid value encountered in sqrt
  return np.sqrt(decode_array(data, **kwargs))
[25]:
# which keywords can be used to access the content?
print(fcontent.keys())
# print the entire content including values of data and
# metadata of the first sweep
# (numpy arrays will not be entirely printed)
print(fcontent["data"][1].keys())
print()
print(fcontent["data"][1]["ingest_data_hdrs"].keys())
print(fcontent["data"][1]["ingest_data_hdrs"]["DB_DBZ"])
print()
print(fcontent["data"][1]["sweep_data"].keys())
print(fcontent["data"][1]["sweep_data"]["DB_DBZ"])
odict_keys(['product_hdr', 'product_type', 'ingest_header', 'nsweeps', 'nrays', 'nbins', 'data_types', 'data', 'raw_product_bhdrs'])
odict_keys(['ingest_data_hdrs', 'record_number', 'first_ray_byte_offset', 'sweep_data', 'ray_offsets'])

odict_keys(['DB_DBZ', 'DB_VEL', 'DB_ZDR', 'DB_KDP', 'DB_PHIDP', 'DB_RHOHV', 'DB_HCLASS'])
OrderedDict([('structure_header', OrderedDict([('structure_identifier', 24), ('format_version', 3), ('bytes_in_structure', 244876), ('flag', 1)])), ('sweep_start_time', datetime.datetime(2013, 11, 25, 10, 55, 3, 541000, tzinfo=datetime.timezone.utc)), ('sweep_number', 1), ('number_rays_per_sweep', 360), ('first_ray_index', 0), ('number_rays_file_expected', 360), ('number_rays_file_written', 360), ('fixed_angle', 0.4998779296875), ('bits_per_bin', 8), ('data_type', 2)])

odict_keys(['DB_DBZ', 'azi_start', 'azi_stop', 'azimuth', 'ele_start', 'ele_stop', 'elevation', 'rbins', 'dtime', 'DB_VEL', 'DB_ZDR', 'DB_KDP', 'DB_PHIDP', 'DB_RHOHV', 'DB_HCLASS'])
[[-20.5 -32.  -32.  ... -32.  -32.  -32. ]
 [-32.  -32.  -32.  ... -32.  -32.  -32. ]
 [-32.  -32.   20.5 ... -32.  -32.  -32. ]
 ...
 [-32.  -32.  -32.  ... -32.  -32.  -32. ]
 [-20.  -32.  -32.  ... -32.  -32.  -32. ]
 [-32.    3.5   6.  ... -32.  -32.  -32. ]]
[26]:
fig = pl.figure(figsize=(10, 10))
swp = fcontent["data"][1]["sweep_data"]
ax, im = wrl.vis.plot_ppi(swp["DB_DBZ"], fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_59_0.png
[27]:
fpath = "sigmet/cor-main131125105503.RAW2049"
f = wrl.util.get_wradlib_data_file(fpath)
fcontent = wrl.io.read_iris(f, keep_old_sweep_data=True)
/home/runner/micromamba-root/envs/wradlib-notebooks/lib/python3.11/site-packages/xradar/io/backends/iris.py:234: RuntimeWarning: invalid value encountered in sqrt
  return np.sqrt(decode_array(data, **kwargs))
[28]:
# which keywords can be used to access the content?
print(fcontent.keys())
# print the entire content including values of data and
# metadata of the first sweep
# (numpy arrays will not be entirely printed)
print(fcontent["data"][1].keys())
print()
print(fcontent["data"][1]["ingest_data_hdrs"].keys())
print(fcontent["data"][1]["ingest_data_hdrs"]["DB_DBZ"])
print()
print(fcontent["data"][1]["sweep_data"].keys())
print(fcontent["data"][1]["sweep_data"]["DB_DBZ"])
odict_keys(['product_hdr', 'product_type', 'ingest_header', 'nsweeps', 'nrays', 'nbins', 'data_types', 'data', 'raw_product_bhdrs'])
odict_keys(['ingest_data_hdrs', 'record_number', 'first_ray_byte_offset', 'sweep_data', 'ray_offsets'])

odict_keys(['DB_DBZ', 'DB_VEL', 'DB_ZDR', 'DB_KDP', 'DB_PHIDP', 'DB_RHOHV', 'DB_HCLASS'])
OrderedDict([('structure_header', OrderedDict([('structure_identifier', 24), ('format_version', 3), ('bytes_in_structure', 244876), ('flag', 1)])), ('sweep_start_time', datetime.datetime(2013, 11, 25, 10, 55, 3, 541000, tzinfo=datetime.timezone.utc)), ('sweep_number', 1), ('number_rays_per_sweep', 360), ('first_ray_index', 0), ('number_rays_file_expected', 360), ('number_rays_file_written', 360), ('fixed_angle', 0.4998779296875), ('bits_per_bin', 8), ('data_type', 2)])

odict_keys(['DB_DBZ', 'DB_VEL', 'DB_ZDR', 'DB_KDP', 'DB_PHIDP', 'DB_RHOHV', 'DB_HCLASS'])
{'data': array([[-20.5, -32. , -32. , ..., -32. , -32. , -32. ],
       [-32. , -32. , -32. , ..., -32. , -32. , -32. ],
       [-32. , -32. ,  20.5, ..., -32. , -32. , -32. ],
       ...,
       [-32. , -32. , -32. , ..., -32. , -32. , -32. ],
       [-20. , -32. , -32. , ..., -32. , -32. , -32. ],
       [-32. ,   3.5,   6. , ..., -32. , -32. , -32. ]]), 'azimuth': array([1.13433838e+00, 2.09014893e+00, 2.99102783e+00, 4.15557861e+00,
       4.94659424e+00, 6.19354248e+00, 7.00378418e+00, 8.10516357e+00,
       9.07196045e+00, 1.00003052e+01, 1.10165405e+01, 1.20629883e+01,
       1.30654907e+01, 1.40652466e+01, 1.50677490e+01, 1.60647583e+01,
       1.70425415e+01, 1.80477905e+01, 1.90805054e+01, 2.00225830e+01,
       2.10690308e+01, 2.20303345e+01, 2.29531860e+01, 2.41369629e+01,
       2.51531982e+01, 2.61529541e+01, 2.69027710e+01, 2.80151367e+01,
       2.90341187e+01, 3.00778198e+01, 3.10693359e+01, 3.20608521e+01,
       3.30139160e+01, 3.40438843e+01, 3.51123047e+01, 3.60598755e+01,
       3.70046997e+01, 3.81555176e+01, 3.91580200e+01, 3.99215698e+01,
       4.10723877e+01, 4.20281982e+01, 4.30664062e+01, 4.41046143e+01,
       4.50494385e+01, 4.59695435e+01, 4.71752930e+01, 4.80789185e+01,
       4.89578247e+01, 5.00180054e+01, 5.10369873e+01, 5.21163940e+01,
       5.29595947e+01, 5.39950562e+01, 5.51046753e+01, 5.59725952e+01,
       5.70574951e+01, 5.79803467e+01, 5.90130615e+01, 6.01171875e+01,
       6.09713745e+01, 6.20452881e+01, 6.30505371e+01, 6.40750122e+01,
       6.50143433e+01, 6.60388184e+01, 6.70056152e+01, 6.82250977e+01,
       6.89666748e+01, 7.00598145e+01, 7.10842896e+01, 7.19082642e+01,
       7.31579590e+01, 7.40148926e+01, 7.50613403e+01, 7.61325073e+01,
       7.70663452e+01, 7.80166626e+01, 7.90521240e+01, 8.00024414e+01,
       8.11285400e+01, 8.20541382e+01, 8.30154419e+01, 8.41113281e+01,
       8.49710083e+01, 8.61108398e+01, 8.70474243e+01, 8.80718994e+01,
       8.90249634e+01, 9.01016235e+01, 9.09832764e+01, 9.20626831e+01,
       9.30761719e+01, 9.39962769e+01, 9.50289917e+01, 9.60644531e+01,
       9.70477295e+01, 9.80502319e+01, 9.90692139e+01, 1.00014038e+02,
       1.01057739e+02, 1.02062988e+02, 1.03035278e+02, 1.04070740e+02,
       1.05070496e+02, 1.06048279e+02, 1.06990356e+02, 1.08053284e+02,
       1.09042053e+02, 1.09984131e+02, 1.11044312e+02, 1.12022095e+02,
       1.13112488e+02, 1.14010620e+02, 1.15057068e+02, 1.16018372e+02,
       1.17026367e+02, 1.18094788e+02, 1.19053345e+02, 1.20061340e+02,
       1.21025391e+02, 1.22091064e+02, 1.22953491e+02, 1.24030151e+02,
       1.25065613e+02, 1.26059875e+02, 1.26990967e+02, 1.28031921e+02,
       1.29053650e+02, 1.30069885e+02, 1.31044922e+02, 1.32061157e+02,
       1.33085632e+02, 1.33989258e+02, 1.35104370e+02, 1.36071167e+02,
       1.37013245e+02, 1.38043213e+02, 1.39045715e+02, 1.40018005e+02,
       1.41045227e+02, 1.42042236e+02, 1.43091431e+02, 1.44047241e+02,
       1.45019531e+02, 1.46079712e+02, 1.47002563e+02, 1.48070984e+02,
       1.49040527e+02, 1.50043030e+02, 1.51056519e+02, 1.52053528e+02,
       1.53056030e+02, 1.54047546e+02, 1.55036316e+02, 1.56041565e+02,
       1.57055054e+02, 1.58041077e+02, 1.59112244e+02, 1.60004883e+02,
       1.61087036e+02, 1.62007141e+02, 1.63094788e+02, 1.64003906e+02,
       1.65044861e+02, 1.66033630e+02, 1.66931763e+02, 1.68104553e+02,
       1.68966980e+02, 1.70109558e+02, 1.71029663e+02, 1.72051392e+02,
       1.73048401e+02, 1.74001465e+02, 1.75036926e+02, 1.76064148e+02,
       1.77157288e+02, 1.77962036e+02, 1.79077148e+02, 1.80016479e+02,
       1.81079407e+02, 1.82043457e+02, 1.83100891e+02, 1.84007263e+02,
       1.84954834e+02, 1.86067200e+02, 1.86970825e+02, 1.88074951e+02,
       1.89044495e+02, 1.90038757e+02, 1.91030273e+02, 1.92082214e+02,
       1.93057251e+02, 1.94051514e+02, 1.95037537e+02, 1.96051025e+02,
       1.97135925e+02, 1.97992859e+02, 1.99143677e+02, 2.00000610e+02,
       2.01003113e+02, 2.02060547e+02, 2.03156433e+02, 2.03919983e+02,
       2.05065308e+02, 2.06023865e+02, 2.07062073e+02, 2.08048096e+02,
       2.09056091e+02, 2.10044861e+02, 2.11105042e+02, 2.11981201e+02,
       2.13096313e+02, 2.13983459e+02, 2.15084839e+02, 2.16109314e+02,
       2.16985474e+02, 2.18191223e+02, 2.18957520e+02, 2.20028687e+02,
       2.21072388e+02, 2.22047424e+02, 2.23038940e+02, 2.24121094e+02,
       2.24953308e+02, 2.26148071e+02, 2.26983032e+02, 2.28076172e+02,
       2.29122620e+02, 2.30067444e+02, 2.31094666e+02, 2.31946106e+02,
       2.33044739e+02, 2.34154358e+02, 2.35066223e+02, 2.36005554e+02,
       2.36997070e+02, 2.38035278e+02, 2.39084473e+02, 2.40070496e+02,
       2.41138916e+02, 2.41924438e+02, 2.43097229e+02, 2.43998108e+02,
       2.45030823e+02, 2.46052551e+02, 2.47060547e+02, 2.47986145e+02,
       2.49082031e+02, 2.50087280e+02, 2.51073303e+02, 2.52009888e+02,
       2.53056335e+02, 2.54113770e+02, 2.55000916e+02, 2.55995178e+02,
       2.57055359e+02, 2.58186951e+02, 2.59087830e+02, 2.60106812e+02,
       2.60936279e+02, 2.62114563e+02, 2.63048401e+02, 2.64103088e+02,
       2.65028687e+02, 2.66022949e+02, 2.67099609e+02, 2.68025208e+02,
       2.69000244e+02, 2.70211487e+02, 2.70977783e+02, 2.72095642e+02,
       2.72927856e+02, 2.74166565e+02, 2.75007019e+02, 2.76072693e+02,
       2.77212524e+02, 2.77866211e+02, 2.79044495e+02, 2.80052490e+02,
       2.81027527e+02, 2.82090454e+02, 2.83021545e+02, 2.83969116e+02,
       2.85032043e+02, 2.85993347e+02, 2.87190857e+02, 2.88110962e+02,
       2.89102478e+02, 2.89942932e+02, 2.91093750e+02, 2.92002869e+02,
       2.93104248e+02, 2.93999634e+02, 2.95073547e+02, 2.96243591e+02,
       2.96976929e+02, 2.98009644e+02, 2.99116516e+02, 3.00036621e+02,
       3.01063843e+02, 3.02077332e+02, 3.02904053e+02, 3.04057617e+02,
       3.04977722e+02, 3.06213684e+02, 3.06968994e+02, 3.08122559e+02,
       3.08935547e+02, 3.10100098e+02, 3.10962524e+02, 3.12132568e+02,
       3.13025208e+02, 3.14066162e+02, 3.15063171e+02, 3.16040955e+02,
       3.17005005e+02, 3.18153076e+02, 3.18996277e+02, 3.19930115e+02,
       3.21045227e+02, 3.21929626e+02, 3.23124390e+02, 3.24069214e+02,
       3.25244751e+02, 3.25923157e+02, 3.27128906e+02, 3.27966614e+02,
       3.29125671e+02, 3.30133667e+02, 3.30993347e+02, 3.32015076e+02,
       3.33050537e+02, 3.34003601e+02, 3.35154419e+02, 3.36000366e+02,
       3.37137451e+02, 3.38167419e+02, 3.38903503e+02, 3.39960938e+02,
       3.41084290e+02, 3.42067566e+02, 3.43050842e+02, 3.44306030e+02,
       3.44921265e+02, 3.46006165e+02, 3.47113037e+02, 3.48057861e+02,
       3.49093323e+02, 3.50038147e+02, 3.50988464e+02, 3.52062378e+02,
       3.52963257e+02, 3.54254150e+02, 3.54962769e+02, 3.56143799e+02,
       3.56934814e+02, 3.58049927e+02, 3.58981018e+02, 2.19726562e-02]), 'rbins': array([664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664, 664,
       664, 664, 664, 664, 664, 664, 664, 664, 664], dtype=int16), 'dtime': array([11, 11, 11, 11, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12,
       12, 12, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13, 13,
       14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 15, 15,
       15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 15, 16, 16, 16, 16,
       16, 16, 16, 16, 16, 16, 16, 16, 16, 16, 16, 17, 17, 17, 17, 17, 17,
       17, 17, 17, 17, 17, 17, 17, 17, 17, 18, 18, 18, 18, 18, 18, 18, 18,
       18, 18, 18, 18, 18, 18, 18, 19, 19, 19, 19, 19, 19, 19, 19, 19, 19,
       19, 19, 19, 19, 19, 20, 20, 20, 20, 20, 20, 20, 20, 20, 20, 20, 20,
       20, 20, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21, 21,
       22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 22, 23,
       23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 23, 24, 24, 24, 24,
       24, 24, 24, 24, 24, 24, 24, 24, 24, 24, 24, 25, 25, 25,  1,  1,  1,
        1,  1,  1,  1,  1,  1,  1,  1,  1,  1,  2,  2,  2,  2,  2,  2,  2,
        2,  2,  2,  2,  2,  2,  2,  2,  3,  3,  3,  3,  3,  3,  3,  3,  3,
        3,  3,  3,  3,  3,  3,  4,  4,  4,  4,  4,  4,  4,  4,  4,  4,  4,
        4,  4,  4,  4,  5,  5,  5,  5,  5,  5,  5,  5,  5,  5,  5,  5,  5,
        5,  5,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,  6,
        7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  7,  8,  8,
        8,  8,  8,  8,  8,  8,  8,  8,  8,  8,  8,  8,  8,  9,  9,  9,  9,
        9,  9,  9,  9,  9,  9,  9,  9,  9,  9,  9, 10, 10, 10, 10, 10, 10,
       10, 10, 10, 10, 10, 10, 10, 10, 10, 11, 11, 11, 11, 11, 11, 11, 11,
       11, 11, 11], dtype=int16), 'azi_start': array([  0.69763184,   1.66442871,   2.45544434,   3.79577637,
         4.37255859,   5.86669922,   6.49291992,   7.71240234,
         8.63525391,   9.48669434,  10.51391602,  11.6015625 ,
        12.61230469,  13.57910156,  14.61730957,  15.61157227,
        16.57287598,  17.578125  ,  18.63830566,  19.52270508,
        20.63781738,  21.55517578,  22.39562988,  23.74145508,
        24.79614258,  25.59265137,  26.27929688,  27.50976562,
        28.54797363,  29.61914062,  30.60791016,  31.61865234,
        32.50305176,  33.58520508,  34.71130371,  35.60668945,
        36.48010254,  37.80944824,  38.60595703,  39.32556152,
        40.63842773,  41.53930664,  42.61047363,  43.68713379,
        44.57702637,  45.4119873 ,  46.8347168 ,  47.64221191,
        48.39477539,  49.52087402,  50.57556152,  51.67419434,
        52.39929199,  53.49243164,  54.7064209 ,  55.44250488,
        56.60705566,  57.42553711,  58.4967041 ,  59.72717285,
        60.43029785,  61.5838623 ,  62.58911133,  63.62731934,
        64.52270508,  65.57189941,  66.5057373 ,  67.8515625 ,
        68.4173584 ,  69.60388184,  70.6640625 ,  71.31774902,
        72.80090332,  73.51501465,  74.54223633,  75.74523926,
        76.63513184,  77.49755859,  78.59069824,  79.47509766,
        80.75500488,  81.59545898,  82.51281738,  83.69934082,
        84.44091797,  85.69885254,  86.58325195,  87.63244629,
        88.53881836,  89.69787598,  90.45593262,  91.62597656,
        92.64770508,  93.47717285,  94.54284668,  95.60852051,
        96.58630371,  97.60253906,  98.61877441,  99.51965332,
       100.60180664, 101.6015625 , 102.56286621, 103.62304688,
       104.63378906, 105.59509277, 106.47399902, 107.60559082,
       108.56140137, 109.46228027, 110.58288574, 111.53869629,
       112.70874023, 113.51623535, 114.59289551, 115.52124023,
       116.54846191, 117.66357422, 118.58642578, 119.61364746,
       120.53100586, 121.66809082, 122.39868164, 123.55773926,
       124.62890625, 125.61218262, 126.48010254, 127.5567627 ,
       128.5949707 , 129.62768555, 130.57250977, 131.61071777,
       132.65991211, 133.47290039, 134.69238281, 135.63720703,
       136.50512695, 137.58178711, 138.5925293 , 139.52636719,
       140.59204102, 141.58630371, 142.61901855, 143.59130859,
       144.53613281, 145.65673828, 146.50268555, 147.62329102,
       148.5736084 , 149.56787109, 150.60058594, 151.59484863,
       152.59460449, 153.58337402, 154.56665039, 155.56640625,
       156.61010742, 157.5604248 , 158.70300293, 159.49401855,
       160.65856934, 161.51550293, 162.68005371, 163.50952148,
       164.58618164, 165.55847168, 166.35498047, 167.70629883,
       168.43139648, 169.70031738, 170.54626465, 171.57348633,
       172.59521484, 173.50158691, 174.56726074, 175.62194824,
       176.67114258, 177.4017334 , 178.64318848, 179.5111084 ,
       180.6427002 , 181.58752441, 182.68615723, 183.515625  ,
       184.41101074, 185.61950684, 186.41052246, 187.63549805,
       188.58032227, 189.5690918 , 190.53588867, 191.6619873 ,
       192.5958252 , 193.57910156, 194.55688477, 195.58410645,
       196.75964355, 197.46276855, 198.78112793, 199.47875977,
       200.49499512, 201.61010742, 202.69226074, 203.33496094,
       204.62585449, 205.53771973, 206.59790039, 207.59216309,
       208.60290527, 209.5916748 , 210.70678711, 211.43737793,
       212.67883301, 213.4588623 , 214.66186523, 215.61218262,
       216.43615723, 217.87536621, 218.39172363, 219.55627441,
       220.62194824, 221.57775879, 222.57202148, 223.73657227,
       224.38476562, 225.7800293 , 226.46118164, 227.62573242,
       228.64196777, 229.63623047, 230.67443848, 231.38305664,
       232.57507324, 233.80554199, 234.63500977, 235.49743652,
       236.49169922, 237.54089355, 238.66149902, 239.6282959 ,
       240.7598877 , 241.33666992, 242.67700195, 243.48999023,
       244.56115723, 245.59936523, 246.62109375, 247.44506836,
       248.66455078, 249.65332031, 250.64208984, 251.5045166 ,
       252.60314941, 253.71826172, 254.48730469, 255.48156738,
       256.60217285, 257.73925781, 258.66760254, 259.70581055,
       260.35949707, 261.71081543, 262.57873535, 263.69384766,
       264.53430176, 265.52307129, 266.67663574, 267.55004883,
       268.47290039, 269.90661621, 270.44494629, 271.6809082 ,
       272.35656738, 273.80126953, 274.4934082 , 275.61950684,
       276.7401123 , 277.22351074, 278.58032227, 279.59655762,
       280.54138184, 281.67297363, 282.53540039, 283.4197998 ,
       284.53491211, 285.46875   , 286.75415039, 287.69348145,
       288.70422363, 289.36889648, 290.6817627 , 291.47827148,
       292.68127441, 293.4942627 , 294.63684082, 295.79040527,
       296.43310547, 297.49328613, 298.72375488, 299.55871582,
       300.62438965, 301.63513184, 302.31079102, 303.59069824,
       304.44213867, 305.85388184, 306.43066406, 307.72705078,
       308.36975098, 309.67712402, 310.41320801, 311.75354004,
       312.5390625 , 313.61572266, 314.60998535, 315.57678223,
       316.50512695, 317.79052734, 318.48815918, 319.35058594,
       320.57006836, 321.36108398, 322.734375  , 323.60778809,
       324.7668457 , 325.33813477, 326.7388916 , 327.41455078,
       328.68896484, 329.75463867, 330.48522949, 331.50146484,
       332.59460449, 333.48449707, 334.74243164, 335.4675293 ,
       336.76391602, 337.62084961, 338.30749512, 339.41162109,
       340.64208984, 341.61437988, 342.59216309, 343.86108398,
       344.33349609, 345.50354004, 346.72851562, 347.53051758,
       348.68408203, 349.55200195, 350.46936035, 351.60644531,
       352.40844727, 353.89709473, 354.41345215, 355.77575684,
       356.35803223, 357.57751465, 358.45092773, 359.54406738]), 'ele_start': array([0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527]), 'elevation': array([0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527]), 'azi_stop': array([  1.57104492,   2.51586914,   3.52661133,   4.51538086,
         5.52062988,   6.52038574,   7.51464844,   8.4979248 ,
         9.50866699,  10.51391602,  11.51916504,  12.52441406,
        13.51867676,  14.5513916 ,  15.51818848,  16.51794434,
        17.51220703,  18.51745605,  19.52270508,  20.52246094,
        21.50024414,  22.50549316,  23.51074219,  24.5324707 ,
        25.51025391,  26.71325684,  27.52624512,  28.52050781,
        29.52026367,  30.53649902,  31.53076172,  32.50305176,
        33.52478027,  34.50256348,  35.51330566,  36.51306152,
        37.52929688,  38.50158691,  39.71008301,  40.51757812,
        41.50634766,  42.51708984,  43.52233887,  44.52209473,
        45.52185059,  46.52709961,  47.51586914,  48.515625  ,
        49.52087402,  50.51513672,  51.49841309,  52.55859375,
        53.51989746,  54.49768066,  55.50292969,  56.50268555,
        57.50793457,  58.53515625,  59.52941895,  60.50720215,
        61.51245117,  62.50671387,  63.51196289,  64.52270508,
        65.50598145,  66.5057373 ,  67.50549316,  68.59863281,
        69.51599121,  70.51574707,  71.5045166 ,  72.4987793 ,
        73.51501465,  74.51477051,  75.58044434,  76.51977539,
        77.49755859,  78.5357666 ,  79.5135498 ,  80.52978516,
        81.5020752 ,  82.51281738,  83.51806641,  84.52331543,
        85.50109863,  86.52282715,  87.51159668,  88.51135254,
        89.5111084 ,  90.50537109,  91.51062012,  92.49938965,
        93.50463867,  94.51538086,  95.51513672,  96.52038574,
        97.50915527,  98.4979248 ,  99.51965332, 100.50842285,
       101.51367188, 102.52441406, 103.50769043, 104.51843262,
       105.50720215, 106.50146484, 107.50671387, 108.50097656,
       109.52270508, 110.50598145, 111.5057373 , 112.50549316,
       113.51623535, 114.50500488, 115.52124023, 116.51550293,
       117.50427246, 118.52600098, 119.52026367, 120.5090332 ,
       121.51977539, 122.51403809, 123.50830078, 124.50256348,
       125.50231934, 126.50756836, 127.50183105, 128.50708008,
       129.5123291 , 130.51208496, 131.51733398, 132.51159668,
       133.51135254, 134.50561523, 135.51635742, 136.50512695,
       137.5213623 , 138.50463867, 139.49890137, 140.50964355,
       141.49841309, 142.49816895, 143.56384277, 144.50317383,
       145.50292969, 146.50268555, 147.50244141, 148.51867676,
       149.50744629, 150.51818848, 151.51245117, 152.51220703,
       153.51745605, 154.51171875, 155.50598145, 156.51672363,
       157.5       , 158.52172852, 159.52148438, 160.51574707,
       161.51550293, 162.4987793 , 163.50952148, 164.49829102,
       165.50354004, 166.50878906, 167.50854492, 168.50280762,
       169.50256348, 170.51879883, 171.51306152, 172.52929688,
       173.50158691, 174.50134277, 175.5065918 , 176.50634766,
       177.64343262, 178.52233887, 179.5111084 , 180.52185059,
       181.51611328, 182.49938965, 183.515625  , 184.49890137,
       185.49865723, 186.51489258, 187.53112793, 188.5144043 ,
       189.50866699, 190.50842285, 191.5246582 , 192.50244141,
       193.51867676, 194.52392578, 195.51818848, 196.51794434,
       197.51220703, 198.52294922, 199.50622559, 200.52246094,
       201.51123047, 202.51098633, 203.62060547, 204.50500488,
       205.50476074, 206.51000977, 207.52624512, 208.50402832,
       209.50927734, 210.49804688, 211.5032959 , 212.52502441,
       213.51379395, 214.50805664, 215.5078125 , 216.60644531,
       217.53479004, 218.50708008, 219.52331543, 220.50109863,
       221.52282715, 222.51708984, 223.50585938, 224.50561523,
       225.52185059, 226.51611328, 227.50488281, 228.52661133,
       229.60327148, 230.49865723, 231.51489258, 232.50915527,
       233.5144043 , 234.50317383, 235.49743652, 236.51367188,
       237.50244141, 238.52966309, 239.50744629, 240.51269531,
       241.51794434, 242.51220703, 243.51745605, 244.50622559,
       245.50048828, 246.5057373 , 247.5       , 248.52722168,
       249.49951172, 250.52124023, 251.5045166 , 252.51525879,
       253.50952148, 254.50927734, 255.51452637, 256.50878906,
       257.50854492, 258.63464355, 259.50805664, 260.5078125 ,
       261.51306152, 262.51831055, 263.51806641, 264.5123291 ,
       265.52307129, 266.52282715, 267.52258301, 268.50036621,
       269.52758789, 270.51635742, 271.51062012, 272.51037598,
       273.49914551, 274.53186035, 275.52062988, 276.52587891,
       277.68493652, 278.50891113, 279.50866699, 280.50842285,
       281.51367188, 282.50793457, 283.50769043, 284.51843262,
       285.5291748 , 286.51794434, 287.62756348, 288.52844238,
       289.50073242, 290.51696777, 291.5057373 , 292.52746582,
       293.52722168, 294.50500488, 295.51025391, 296.69677734,
       297.52075195, 298.52600098, 299.50927734, 300.51452637,
       301.5032959 , 302.51953125, 303.49731445, 304.52453613,
       305.51330566, 306.57348633, 307.50732422, 308.51806641,
       309.50134277, 310.52307129, 311.51184082, 312.51159668,
       313.51135254, 314.51660156, 315.51635742, 316.50512695,
       317.50488281, 318.515625  , 319.50439453, 320.50964355,
       321.52038574, 322.49816895, 323.5144043 , 324.53063965,
       325.72265625, 326.50817871, 327.5189209 , 328.51867676,
       329.56237793, 330.51269531, 331.50146484, 332.52868652,
       333.50646973, 334.52270508, 335.56640625, 336.53320312,
       337.51098633, 338.71398926, 339.49951172, 340.51025391,
       341.52648926, 342.52075195, 343.50952148, 344.75097656,
       345.5090332 , 346.50878906, 347.49755859, 348.58520508,
       349.50256348, 350.52429199, 351.50756836, 352.51831055,
       353.51806641, 354.61120605, 355.51208496, 356.51184082,
       357.51159668, 358.52233887, 359.5111084 ,   0.49987793]), 'ele_stop': array([0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527,
       0.47790527, 0.47790527, 0.47790527, 0.47790527, 0.47790527])}
[29]:
fig = pl.figure(figsize=(10, 10))
swp = fcontent["data"][1]["sweep_data"]
ax, im = wrl.vis.plot_ppi(swp["DB_DBZ"]["data"], fig=fig, proj="cg")
../../_images/notebooks_fileio_wradlib_radar_formats_62_0.png

Furuno SCN and SCNX

Please refer to Furuno Operator Manual for a description of the file format.

Note

Since \(\omega radlib\) version 1.15 a Furuno reader wradlib.io.furuno.open_furuno_dataset() based on Xarray is available. Please read the more indepth notebook wradlib_furuno_backend.

From \(\omega radlib\) version 1.19 Furuno reading code is imported from xradar-package whenever and wherever necessary.

The use of the above xarray backend is recommended , but for closer examination of the header structure the data can be opened using:

[30]:
fpath = "furuno/0080_20210730_160000_01_02.scn.gz"
f = wrl.util.get_wradlib_data_file(fpath)
with wrl.io.furuno.FurunoFile(f, loaddata=True) as scn_file:
    display(scn_file.header)
/tmp/ipykernel_4362/918775053.py:3: FutureWarning: FurunoFile class has been moved to xradar-package. Importing from wradlib will be removed in v2.0.
  with wrl.io.furuno.FurunoFile(f, loaddata=True) as scn_file:
OrderedDict([('size_of_header', 80),
             ('format_version', 3),
             ('dpu_log_time', datetime.datetime(2021, 7, 30, 18, 0)),
             ('latitude', 47.07734000000001),
             ('longitude', 15.44729),
             ('altitude', 407.9),
             ('antenna_rotation_speed', 40),
             ('prf_1', 10000),
             ('prf_2', 8000),
             ('noise_level_pulse_modulation_h', -4983),
             ('noise_level_frequency_modulation_h', -7116),
             ('number_sweep_direction_data', 1376),
             ('number_range_direction_data', 602),
             ('resolution_range_direction', 5000),
             ('radar_constant_h', 9.575470000000001e-17),
             ('radar_constant_v', 9.144499999999999e-17),
             ('azimuth_offset', 15595),
             ('scan_start_time', datetime.datetime(2021, 7, 30, 16, 0)),
             ('record_item', 511),
             ('tx_pulse_blind_length', 163),
             ('tx_pulse_specification', 4)])
[31]:
fig = pl.figure(figsize=(10, 10))
dbzh = scn_file.data["DBZH"]
az = scn_file.data["azimuth"] / 100
rres = scn_file.header["resolution_range_direction"] / 100
r = np.arange(rres / 2, dbzh.shape[-1] * rres, rres)
ax, im = wrl.vis.plot_ppi(dbzh, az=az, r=r, fig=fig, proj="cg")
/home/runner/micromamba-root/envs/wradlib-notebooks/lib/python3.11/site-packages/xarray/plot/dataarray_plot.py:2313: UserWarning: The input coordinates to pcolormesh are interpreted as cell centers, but are not monotonically increasing or decreasing. This may lead to incorrectly calculated cell edges, in which case, please supply explicit cell edges to pcolormesh.
  primitive = ax.pcolormesh(x, y, z, **kwargs)
../../_images/notebooks_fileio_wradlib_radar_formats_65_1.png
[32]:
fpath = "furuno/2006_20220324_000000_000.scnx.gz"
f = wrl.util.get_wradlib_data_file(fpath)
with wrl.io.furuno.FurunoFile(f, loaddata=True) as scnx_file:
    display(scnx_file.header)
/tmp/ipykernel_4362/1871657942.py:3: FutureWarning: FurunoFile class has been moved to xradar-package. Importing from wradlib will be removed in v2.0.
  with wrl.io.furuno.FurunoFile(f, loaddata=True) as scnx_file:
OrderedDict([('size_of_header', 156),
             ('format_version', 10),
             ('scan_start_time', datetime.datetime(2022, 3, 24, 0, 0, 1)),
             ('scan_stop_time', datetime.datetime(2022, 3, 24, 0, 0, 29)),
             ('time_zone', 100),
             ('product_number', 2006),
             ('model_type', 4),
             ('latitude', 5355478),
             ('longitude', 1324397),
             ('altitude', 3800),
             ('azimuth_offset', 27160),
             ('tx_frequency', 943250),
             ('polarization_mode', 2),
             ('antenna_gain_h', 337),
             ('antenna_gain_v', 338),
             ('half_power_beam_width_h', 275),
             ('half_power_beam_width_v', 276),
             ('tx_power_h', 724),
             ('tx_power_v', 708),
             ('radar_constant_h', -1294),
             ('radar_constant_v', -1296),
             ('noise_power_short_pulse_h', -604),
             ('noise_power_long_pulse_h', -767),
             ('threshold_power_short_pulse', -564),
             ('threshold_power_long_pulse', -785),
             ('tx_pulse_specification', 4),
             ('prf_mode', 2),
             ('prf_1', 13000),
             ('prf_2', 10400),
             ('prf_3', 0),
             ('nyquist_velocity', 413),
             ('sample_number', 165),
             ('tx_pulse_blind_length', 263),
             ('short_pulse_width', 100),
             ('short_pulse_modulation_bandwidth', 0),
             ('long_pulse_width', 5000),
             ('long_pulse_modulation_bandwidth', 200),
             ('pulse_switch_point', 7650),
             ('observation_mode', 3),
             ('antenna_rotation_speed', 21),
             ('number_sweep_direction_data', 722),
             ('number_range_direction_data', 936),
             ('resolution_range_direction', 75),
             ('current_scan_number', 0),
             ('total_number_scans_volume', 9),
             ('rainfall_intensity_estimation_method', 4),
             ('z_r_coefficient_b', 2000),
             ('z_r_coefficient_beta', 160),
             ('kdp_r_coefficient_a', 1960),
             ('kdp_r_coefficient_b', 815),
             ('kdp_r_coefficient_c', 120),
             ('zh_attenuation_correction_method', 3),
             ('zh_attenuation_coefficient_b1', 293),
             ('zh_attenuation_coefficient_b2', 1101),
             ('zh_attenuation_coefficient_d1', 298),
             ('zh_attenuation_coefficient_d2', 1293),
             ('air_attenuation_one_way', 10),
             ('output_threshold_rain', 0),
             ('record_item', 33279),
             ('signal_processing_flag', 3),
             ('clutter_reference_file',
              datetime.datetime(2021, 7, 19, 12, 0))])
[33]:
fig = pl.figure(figsize=(10, 10))
dbzh = scnx_file.data["DBZH"]
az = scnx_file.data["azimuth"] / 100
rres = scnx_file.header["resolution_range_direction"] / 100
r = np.arange(rres / 2, dbzh.shape[-1] * rres, rres)
ax, im = wrl.vis.plot_ppi(dbzh, az=az, r=r, fig=fig, proj="cg")
/home/runner/micromamba-root/envs/wradlib-notebooks/lib/python3.11/site-packages/xarray/plot/dataarray_plot.py:2313: UserWarning: The input coordinates to pcolormesh are interpreted as cell centers, but are not monotonically increasing or decreasing. This may lead to incorrectly calculated cell edges, in which case, please supply explicit cell edges to pcolormesh.
  primitive = ax.pcolormesh(x, y, z, **kwargs)
../../_images/notebooks_fileio_wradlib_radar_formats_67_1.png

OPERA BUFR

WARNING \(\omega radlib\) does currently not support the BUFR format!

The Binary Universal Form for the Representation of meteorological data (BUFR) is a binary data format maintained by the World Meteorological Organization (WMO).

The BUFR format was adopted by OPERA for the representation of weather radar data. A BUFR file consists of a set of descriptors which contain all the relevant metadata and a data section. The descriptors are identified as a tuple of three integers. The meaning of these tupels is described in the so-called BUFR tables. There are generic BUFR tables provided by the WMO, but it is also possible to define so called local tables - which was done by the OPERA consortium for the purpose of radar data representation.

If you want to use BUFR files together with \(\omega radlib\), we recommend that you check out the OPERA webpage where you will find software for BUFR decoding. In particular, you might want to check out this tool which seems to support the conversion of OPERA BUFR files to ODIM_H5 (which is supported by \(\omega radlib\)). However, you have to build it yourself.

It would be great if someone could add a tutorial on how to use OPERA BUFR software together with \(\omega radlib\)!