{ "cells": [ { "cell_type": "markdown", "id": "59447ad6-ac47-494e-b696-4335b36b205b", "metadata": {}, "source": [ "# CfRadial1 to CfRadial2 - A data model transformation\n", "\n", "In this notebook we show how to transform the CfRadial1 Data model to a CfRadial2 representation.\n", "\n", "We use some internal functions to show how xradar is working inside.\n", "\n", "Within this notebook we reference to the [CfRadial2.1 draft](https://github.com/NCAR/CfRadial/tree/master/docs). As long as the FM301 WMO standard is not finalized we will rely on the drafts presented." ] }, { "cell_type": "code", "execution_count": null, "id": "6f96b5d8-2b96-4fd7-b8ba-166c34a8dcd2", "metadata": {}, "outputs": [], "source": [ "import os\n", "import xarray as xr\n", "import xradar as xd\n", "import datatree as xt\n", "from open_radar_data import DATASETS" ] }, { "cell_type": "markdown", "id": "33d50be4-dfe5-4d99-a936-67a9a76bac94", "metadata": {}, "source": [ "## Download\n", "\n", "Fetching CfRadial1 radar data file from [open-radar-data](https://github.com/openradar/open-radar-data) repository." ] }, { "cell_type": "code", "execution_count": null, "id": "d3c6d408-5ab2-43c3-afd1-b3a703ef3b24", "metadata": {}, "outputs": [], "source": [ "filename = DATASETS.fetch(\"cfrad.20080604_002217_000_SPOL_v36_SUR.nc\")" ] }, { "cell_type": "markdown", "id": "b987dcfd-5105-4483-932e-71b8002e5f09", "metadata": {}, "source": [ "## Open CfRadial1 file using xr.open_dataset\n", "\n", "Making use of the xarray `netcdf4` backend. We get back all data and metadata in one single CfRadial1 Dataset." ] }, { "cell_type": "code", "execution_count": null, "id": "0e0a7a18-3a4c-4940-96e0-6059f0b3da48", "metadata": {}, "outputs": [], "source": [ "ds = xr.open_dataset(filename, engine=\"netcdf4\")\n", "with xr.set_options(\n", " display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000\n", "):\n", " display(ds.load())" ] }, { "cell_type": "markdown", "id": "b37bc43e", "metadata": {}, "source": [ "## Extract CfRadial2 Groups and Subgroups\n", "\n", "Now as we have the CfRadial1 Dataset we can work towards extracting the CfRadial2 groups and subgroups." ] }, { "cell_type": "markdown", "id": "e1d1b525", "metadata": {}, "source": [ "### Extract CfRadial2 Root-Group\n", "\n", "The following sections present the details of the information in the top-level (root) group of the\n", "data set.\n", "\n", "We use a convenience function to extract the CfRadial2 root group from the CfRadial1 Dataset. We can call this function with one kwarg:\n", "\n", "- `optional=False` - only mandatory data and metadata is imported, defaults to True" ] }, { "cell_type": "markdown", "id": "f4f3ff2f", "metadata": {}, "source": [ "#### optional=True" ] }, { "cell_type": "code", "execution_count": null, "id": "e4b1557d", "metadata": {}, "outputs": [], "source": [ "root = xd.io.backends.cfradial1._get_required_root_dataset(ds)\n", "with xr.set_options(\n", " display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000\n", "):\n", " display(root.load())" ] }, { "cell_type": "markdown", "id": "7459f970", "metadata": {}, "source": [ "#### optional=False" ] }, { "cell_type": "code", "execution_count": null, "id": "fff2a34d", "metadata": {}, "outputs": [], "source": [ "root = xd.io.backends.cfradial1._get_required_root_dataset(ds, optional=False)\n", "with xr.set_options(\n", " display_expand_data_vars=True, display_expand_attrs=True, display_max_rows=1000\n", "):\n", " display(root)" ] }, { "cell_type": "markdown", "id": "65cf9049", "metadata": {}, "source": [ "### Extract Root-Group metadata groups\n", "\n", "The Cfradial2 Data Model has a notion of root group metadata groups. Those groups provide additional metadata covering other aspects of the radar system.\n", "\n", "#### The radar_parameters sub-group\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": null, "id": "b992ec33", "metadata": {}, "outputs": [], "source": [ "display(xd.model.radar_parameters_subgroup)" ] }, { "cell_type": "markdown", "id": "fd5589ce", "metadata": {}, "source": [ "Again we use a convenience function to extract the group." ] }, { "cell_type": "code", "execution_count": null, "id": "04af8243", "metadata": {}, "outputs": [], "source": [ "radar_parameters = xd.io.backends.cfradial1._get_subgroup(\n", " ds, xd.model.radar_parameters_subgroup\n", ")\n", "display(radar_parameters.load())" ] }, { "cell_type": "markdown", "id": "c926508e", "metadata": {}, "source": [ "#### The radar_calibration sub-group\n", "\n", "For a radar, a different calibration is required for each pulse width. Therefore the calibration\n", "variables are arrays. If only one calibration is available it is squeezed by the reader." ] }, { "cell_type": "code", "execution_count": null, "id": "8cb5f197", "metadata": {}, "outputs": [], "source": [ "display(xd.model.radar_calibration_subgroup)" ] }, { "cell_type": "markdown", "id": "6bcaa501", "metadata": {}, "source": [ "Again we use a convenience function to extract the group." ] }, { "cell_type": "code", "execution_count": null, "id": "9f61ce7f", "metadata": {}, "outputs": [], "source": [ "radar_calibration = xd.io.backends.cfradial1._get_radar_calibration(ds)\n", "with xr.set_options(display_expand_data_vars=True):\n", " if radar_calibration:\n", " display(radar_calibration.load())" ] }, { "cell_type": "markdown", "id": "2ca0320a", "metadata": {}, "source": [ "#### The georeference_correction sub-group\n", "\n", "The following additional variables are used to quantify errors in the georeference data for moving\n", "platforms. These are constant for a volume." ] }, { "cell_type": "code", "execution_count": null, "id": "9e76d49e", "metadata": {}, "outputs": [], "source": [ "display(xd.model.georeferencing_correction_subgroup)" ] }, { "cell_type": "markdown", "id": "c7739019", "metadata": {}, "source": [ "Again we use a convenience function to extract the group." ] }, { "cell_type": "code", "execution_count": null, "id": "a42c442f", "metadata": {}, "outputs": [], "source": [ "georeference_correction = xd.io.backends.cfradial1._get_subgroup(\n", " ds, xd.model.georeferencing_correction_subgroup\n", ")\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(georeference_correction.load())" ] }, { "cell_type": "markdown", "id": "042da0a9", "metadata": {}, "source": [ "### Sweep groups\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": null, "id": "70224a48", "metadata": {}, "outputs": [], "source": [ "root.sweep_group_name" ] }, { "cell_type": "markdown", "id": "7f4f2c1c-5961-43f6-8c40-cbe796dd12dd", "metadata": {}, "source": [ "Again we use a convenience function to extract the different sweep groups. We can call this function with kwargs:\n", "\n", "- `optional=False` - only mandatory data and metadata is imported, defaults to `True`\n", "- `first_dim=\"time` - return first dimension as `time`, defaults to`auto` (return either as `azimuth` (PPI) or `elevation` (RHI)to `time`\n", "- `site_coords=False` - do not add radar site coordinates to the Sweep-Dataset, defaults to `True`\n", "\n", "#### Examining first sweep with default kwargs." ] }, { "cell_type": "code", "execution_count": null, "id": "cf93cd38-0f00-4dec-b7af-4e0ca6a111d1", "metadata": {}, "outputs": [], "source": [ "sweeps = xd.io.backends.cfradial1._get_sweep_groups(ds)\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(sweeps[\"sweep_0\"])" ] }, { "cell_type": "markdown", "id": "38b326e1", "metadata": {}, "source": [ "#### Examining first sweep with `optional=False`" ] }, { "cell_type": "code", "execution_count": null, "id": "729e1ee3", "metadata": {}, "outputs": [], "source": [ "sweeps = xd.io.backends.cfradial1._get_sweep_groups(ds, optional=False)\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(sweeps[\"sweep_0\"])" ] }, { "cell_type": "markdown", "id": "32652409", "metadata": {}, "source": [ "#### `optional=False` and `site_coords=False`" ] }, { "cell_type": "code", "execution_count": null, "id": "deaaa18e", "metadata": {}, "outputs": [], "source": [ "sweeps = xd.io.backends.cfradial1._get_sweep_groups(\n", " ds, optional=False, site_coords=False\n", ")\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(sweeps[\"sweep_0\"])" ] }, { "cell_type": "markdown", "id": "65a43c95", "metadata": {}, "source": [ "#### `optional=False`, `site_coords=True` and `first_dim=\"auto\"`" ] }, { "cell_type": "code", "execution_count": null, "id": "ef477a6d", "metadata": {}, "outputs": [], "source": [ "sweeps = xd.io.backends.cfradial1._get_sweep_groups(\n", " ds, optional=False, site_coords=False, first_dim=\"time\"\n", ")\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(sweeps[\"sweep_0\"])" ] }, { "cell_type": "markdown", "id": "e320f119", "metadata": {}, "source": [ "## Read as CfRadial2 data representation\n", "\n", "xradar provides two easy ways to retrieve the CfRadial1 data as CfRadial2 groups.\n", "\n", "### DataTree\n", "\n", "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." ] }, { "cell_type": "code", "execution_count": null, "id": "efc02939", "metadata": {}, "outputs": [], "source": [ "dtree = xd.io.open_cfradial1_datatree(filename)\n", "with xr.set_options(display_expand_data_vars=True, display_expand_attrs=True):\n", " display(dtree)" ] }, { "cell_type": "markdown", "id": "666b41bd", "metadata": {}, "source": [ "Each DataTree-node itself represents another DataTree." ] }, { "cell_type": "code", "execution_count": null, "id": "426b4256", "metadata": {}, "outputs": [], "source": [ "display(dtree[\"radar_parameters\"].load())" ] }, { "cell_type": "code", "execution_count": null, "id": "04083ab1", "metadata": {}, "outputs": [], "source": [ "with xr.set_options(display_expand_data_vars=True):\n", " display(dtree[\"sweep_7\"].load())" ] }, { "cell_type": "markdown", "id": "97ee3558-d138-4669-8bcf-91cb493669af", "metadata": {}, "source": [ "#### Roundtrip with `to_netcdf`" ] }, { "cell_type": "markdown", "id": "55515cf9", "metadata": {}, "source": [ "Write DataTree to netCDF4 file, reopen and compare with source. This just tets if roundtripping the DataTree works." ] }, { "cell_type": "code", "execution_count": null, "id": "65cd9185", "metadata": {}, "outputs": [], "source": [ "outfile = \"test_dtree.nc\"\n", "if os.path.exists(outfile):\n", " os.unlink(outfile)\n", "dtree.to_netcdf(outfile)" ] }, { "cell_type": "code", "execution_count": null, "id": "5921f261", "metadata": {}, "outputs": [], "source": [ "dtree2 = xt.open_datatree(outfile)\n", "with xr.set_options(display_expand_data_vars=True, display_expand_attrs=True):\n", " display(dtree2)" ] }, { "cell_type": "code", "execution_count": null, "id": "bf23b4ce", "metadata": {}, "outputs": [], "source": [ "for grp in dtree.groups:\n", " print(grp)\n", " xr.testing.assert_equal(dtree[grp].ds, dtree2[grp].ds)" ] }, { "cell_type": "markdown", "id": "54aaa127", "metadata": {}, "source": [ "#### Roundtrip with `xradar.io.to_cfradial2`" ] }, { "cell_type": "code", "execution_count": null, "id": "ca7ba21f", "metadata": {}, "outputs": [], "source": [ "dtree3 = xd.io.open_cfradial1_datatree(filename)" ] }, { "cell_type": "code", "execution_count": null, "id": "003efee0", "metadata": {}, "outputs": [], "source": [ "display(dtree3)" ] }, { "cell_type": "code", "execution_count": null, "id": "42b00295-08c6-4b03-b293-69fb141239f7", "metadata": {}, "outputs": [], "source": [ "outfile = \"test_cfradial2.nc\"\n", "if os.path.exists(outfile):\n", " os.unlink(outfile)\n", "xd.io.to_cfradial2(dtree3, outfile)" ] }, { "cell_type": "code", "execution_count": null, "id": "de63c0fc-7d9f-443d-b802-e1c903930142", "metadata": {}, "outputs": [], "source": [ "dtree4 = xt.open_datatree(\"test_cfradial2.nc\")\n", "with xr.set_options(display_expand_data_vars=True, display_expand_attrs=True):\n", " display(dtree4)" ] }, { "cell_type": "code", "execution_count": null, "id": "8753ef39-f4eb-48c8-b084-e3e50aee337a", "metadata": {}, "outputs": [], "source": [ "for grp in dtree3.groups:\n", " print(grp)\n", " xr.testing.assert_equal(dtree3[grp].ds, dtree4[grp].ds)" ] }, { "cell_type": "markdown", "id": "fa2c9b84-ce4a-46fb-bd2d-54bcc683f4e1", "metadata": {}, "source": [ "### Datasets\n", "\n", "Using xarray.open_dataset and the cfradial1-backend we can easily load specific groups side-stepping the DataTree. Can be parameterized using kwargs." ] }, { "cell_type": "code", "execution_count": null, "id": "e5e50281-6f8b-4400-b4c9-fb17d07eac24", "metadata": {}, "outputs": [], "source": [ "ds = xr.open_dataset(filename, group=\"sweep_1\", engine=\"cfradial1\", first_dim=\"time\")\n", "with xr.set_options(display_expand_data_vars=True):\n", " display(ds.load())" ] }, { "cell_type": "code", "execution_count": null, "id": "679d0844-6852-43c4-a2ee-78ccd7d8e943", "metadata": {}, "outputs": [], "source": [ "ds = xr.open_dataset(filename, group=\"radar_parameters\", engine=\"cfradial1\")\n", "display(ds.load())" ] }, { "cell_type": "markdown", "id": "84b2975f", "metadata": {}, "source": [ "## Conclusion\n", "\n", "CfRadial1 and CfRadial2 are based on the same principles with slightly different data representation. Nevertheless the conversion is relatively straighforward as has been shown here.\n", "\n", "As the implementation with the cfradial1 xarray backend on one hand and the DataTree on the other hand is very versatile users can pick the most usable approach for their workflows.\n" ] } ], "metadata": { "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.10.6" } }, "nbformat": 4, "nbformat_minor": 5 }