Multipanel hvPlot

Interactive Radar Visualization

Overview

Within this cookbook, we will detail how to create interactive plots of radar data!

  1. Reading data with Xradar

  2. Creating your first interactive figure with Xradar + hvPlot

  3. Combining your plots into a single dashboard

  4. Filtering and Checking Data Quality

  5. Create a Dashboard to Analyze ZDR Bias

Imports

import xradar as xd
import fsspec
import pyart
from open_radar_data import DATASETS
import hvplot.xarray
import holoviews as hv
import panel as pn

hv.extension("bokeh")
## You are using the Python ARM Radar Toolkit (Py-ART), an open source
## library for working with weather radar data. Py-ART is partly
## supported by the U.S. Department of Energy as part of the Atmospheric
## Radiation Measurement (ARM) Climate Research Facility, an Office of
## Science user facility.
##
## If you use this software to prepare a publication, please cite:
##
##     JJ Helmus and SM Collis, JORS 2016, doi: 10.5334/jors.119

Prerequisites

It is recommended that you are familiar with working with weather radar data, the core data structures, and the basics of reading in different radar datasets.

Concepts

Importance

Notes

Intro to Cartopy

Necessary

Xradar User Guide: Plot a PPI

Necessary

Understanding of NetCDF

Helpful

Familiarity with metadata structure

  • Time to learn: 30 minutes

Reading Data with Xradar

While we have focused much of the content around the Python ARM Radar Toolkit (Py-ART), Xradar is another helpful package we can use to work with this in xarray!

Here, we use data from the Colombian weather radar network, using some remote access tools such as fsspec.

fs = fsspec.filesystem("s3", anon=True)
files = sorted(fs.glob("s3-radaresideam/l2_data/2022/08/09/Carimagua/CAR22080919*"))
files
['s3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809190003.RAWDSVV',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809190315.RAWDSW0',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809190401.RAWDSW3',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809190505.RAWDSW8',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191003.RAWDSWM',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191314.RAWDSWT',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191400.RAWDSWX',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191504.RAWDSX2',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191817.RAWDSX7',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809191903.RAWDSXA',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192005.RAWDSXF',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192316.RAWDSXM',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192402.RAWDSXR',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192504.RAWDSXW',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192815.RAWDSY1',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809192901.RAWDSY5',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809193004.RAWDSYA',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809193315.RAWDSYG',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809193402.RAWDSYK',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809193816.RAWDSZ0',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809193902.RAWDSZ4',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194003.RAWDSZ7',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194315.RAWDSZA',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194401.RAWDSZE',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194504.RAWDSZJ',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194816.RAWDSZP',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809194902.RAWDSZU',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195003.RAWDSZY',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195315.RAWDT04',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195401.RAWDT07',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195504.RAWDT0D',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195815.RAWDT0J',
 's3-radaresideam/l2_data/2022/08/09/Carimagua/CAR220809195900.RAWDT0M']

Once we have our files of interest, we can load one in. Let’s read in a single file - for example, we are interested in the 8th file (index=7)

task_file = fsspec.open_local(f'simplecache::s3://{files[7]}',
                                s3={'anon': True},
                                filecache={'cache_storage': '.'})

Read the data in xradar

We use xradar, an open-source toolkit to read weather radar data and load into the Xarray data structure. The data format here is an IRIS file, so we use the open_iris_datatree reader.

radar = xd.io.open_iris_datatree(task_file).xradar.georeference()
radar
<xarray.DatasetView> Size: 248B
Dimensions:              (sweep: 1)
Dimensions without coordinates: sweep
Data variables:
    volume_number        int64 8B 0
    platform_type        <U5 20B 'fixed'
    instrument_type      <U5 20B 'radar'
    time_coverage_start  <U20 80B '2022-08-09T19:15:05Z'
    time_coverage_end    <U20 80B '2022-08-09T19:16:21Z'
    longitude            float64 8B -71.33
    altitude             float64 8B 206.0
    latitude             float64 8B 4.564
    sweep_group_name     (sweep) int64 8B 0
    sweep_fixed_angle    (sweep) float64 8B 0.5
Attributes:
    Conventions:      None
    instrument_name:  carimagua-radar
    version:          None
    title:            None
    institution:      None
    references:       None
    source:           Sigmet
    history:          None
    comment:          Primera tarea del modo precipitacion / 0.5
    scan_name:        SURVP       

Creating Your First Interactive Figure with Xradar and hvPlot

hvPlot is helpful tool when working with interactive visualizions! It is a tool built on top of several other packages, that we can use with Xarray.

By default, this visualization plots azimuth along the y-axis and range along the y-axis. While this is desired in certain cases, we cannot gather much spatial information from this.

ref = radar["sweep_0"].DBZH.hvplot.quadmesh(cmap='pyart_ChaseSpectral',
                                            title='Horizontal Reflectivity (dBZ)',
                                            clim=(-20,60))
ref

Refining Our Plot - Recreating a Plan Position Indicator (PPI)

We instead would like to create a Plan Position Indicator (PPI) plot. Since we already georeferenced the dataset, we set x/y to be x and y, or the distance away from the radar, as well as tuning some additional parameters. We set rasterize=True to lazily load in the data, which renders the plot more quickly and increases resolution as we zoom in.

ref = radar["sweep_0"].DBZH.hvplot.quadmesh(x='x',
                                           y='y',
                                           cmap='pyart_ChaseSpectral',
                                           clabel='Horizontal Reflectivity (dBZ)',
                                           title=f'Horizontal Reflectivity \n {radar.attrs["instrument_name"]} Radar',
                                           clim=(-20, 60),
                                           height=500,
                                             rasterize=True,
                                           width=600,)

ref
zdr = radar["sweep_0"].ZDR.hvplot.quadmesh(x='x',
                                             y='y',
                                             cmap='pyart_ChaseSpectral',
                                             clabel='Differential Reflectivity (dB)',
                                             title=f'Horizontal Reflectivity \n {radar.attrs["instrument_name"]} Radar',
                                             clim=(-1, 6),
                                             height=500,
                                             rasterize=True,
                                             width=600,)
zdr

Combining your plots into a single dashboard

You can combine plots using the + syntax to add plots side-by-side, or * to add them to the same plot. For example, let’s combine our reflectivity and velocity plot.

(ref + zdr).cols(1)

Filtering and Checking Data Quality

We can also filter our data - notice the low values in both differential reflectivity and reflectivity. We can mask these out using Xarray!

# Subset our first sweep
ds = radar["sweep_0"].to_dataset()

Determine Mask Thresholds

Let’s determine our thresholds for filtering the data, using histograms! These are available through hvPlot, using the .hvplot.hist() extension.

zdr_hist = ds.ZDR.hvplot.hist()
ref_hist = ds.DBZH.hvplot.hist()
(zdr_hist + ref_hist).cols(1)

Notice how we have very low values for both fields, which we can threshold using:

  • Differential Reflectivity < -5

  • Horizontal Reflectivity < -32

ds = ds.where((ds.ZDR >= -5) & (ds.DBZH != -32))
ds
<xarray.Dataset> Size: 112MB
Dimensions:            (azimuth: 720, range: 994)
Coordinates:
    elevation          (azimuth) float64 6kB 0.4779 0.4779 ... 0.4779 0.4779
    time               (azimuth) datetime64[ns] 6kB 2022-08-09T19:15:56.89100...
  * range              (range) float32 4kB 1e+03 1.3e+03 ... 2.986e+05 2.989e+05
    longitude          float64 8B -71.33
    latitude           float64 8B 4.564
    altitude           float64 8B 206.0
  * azimuth            (azimuth) float32 3kB 0.03571 0.5795 ... 359.0 359.6
    crs_wkt            int64 8B 0
    x                  (azimuth, range) float64 6MB 0.6231 0.8101 ... -2.191e+03
    y                  (azimuth, range) float64 6MB 999.9 1.3e+03 ... 2.987e+05
    z                  (azimuth, range) float64 6MB 214.4 216.9 ... 7.948e+03
Data variables: (12/17)
    DBTH               (azimuth, range) float64 6MB 28.5 34.0 44.0 ... nan nan
    DBZH               (azimuth, range) float64 6MB 26.0 27.5 28.5 ... nan nan
    VRADH              (azimuth, range) float64 6MB 9.863 9.863 ... nan nan
    WRADH              (azimuth, range) float64 6MB 0.5205 0.5465 ... nan nan
    ZDR                (azimuth, range) float64 6MB 1.125 1.125 ... nan nan
    KDP                (azimuth, range) float64 6MB 0.04458 0.04458 ... nan nan
    ...                 ...
    DB_DBZE8           (azimuth, range) float64 6MB 27.5 29.0 30.0 ... nan nan
    sweep_mode         (azimuth, range) object 6MB 'azimuth_surveillance' ......
    sweep_number       (azimuth, range) float64 6MB 0.0 0.0 0.0 ... nan nan nan
    prt_mode           (azimuth, range) object 6MB 'not_set' 'not_set' ... nan
    follow_mode        (azimuth, range) object 6MB 'not_set' 'not_set' ... nan
    sweep_fixed_angle  (azimuth, range) float64 6MB 0.5 0.5 0.5 ... nan nan nan

Double Check our Filtered Data

Let’s double check that our filtering worked - notice the new, more representative distributions!

zdr_hist = ds.ZDR.hvplot.hist()
ref_hist = ds.DBZH.hvplot.hist()
(zdr_hist + ref_hist).cols(1)

Create a Dashboard to Analyze ZDR Bias

A common data quality check is differential reflectivity bias. This value should be around 0 for low values of horizontal reflectivity. We use a few steps here to create this visualization

  • Unstack the dataset so we are left with a single dimension - the single range gate (single points)

  • Create histograms (.hist) and a 2-dimensional histogram (.hexbin) to visualize the data

  • Stack these into single view using gridspec

ds = ds.stack({"gate": {"azimuth", "range"}}).reset_index("gate")

hist_dbz = ds.hvplot.hist("DBZH",
                          width=500,
                          height=200,
                          title="Horizontal Reflectivity Distribution",)
hist_zdr = ds.hvplot.hist("ZDR",
                          height=400,
                          title="Differential Reflectivity Distribution",
                         ).opts(invert_axes=True)
hexbin = ds.hvplot.hexbin(x='DBZH',
                          y='ZDR',
                          title='Reflectivity vs. Differential Reflectivity Distribution',
                          width=500,
                          height=400) *  hv.HLine(0,
                                                  label='Differential Reflectivity = 0 Line').opts(color='red',
                                                           line_width=1)

gspec = pn.GridSpec(width=800, height=400)

gspec[0,   0:2  ] = hist_dbz
gspec[1:3,   0:2  ] = hexbin
gspec[1:3,   2  ] = hist_zdr

gspec

Summary

Within this notebook, we covered how to use interactive visualizations with your weather radar data, including applications to checking data quality.

What’s Next?

Next, we will continue to explore methods of cleaning and visualizing data!

Resources and References