Py-ART Corrections
Overview
Within this notebook, we will cover:
Intro to radar aliasing.
Calculation of velocity texture using Py-ART
Dealiasing the velocity field
Prerequisites
Concepts |
Importance |
Notes |
---|---|---|
Helpful |
Basic features |
|
Helpful |
Basic plotting |
|
Helpful |
Basic arrays |
Time to learn: 30 minutes
Imports
import os
import warnings
import cartopy.crs as ccrs
import matplotlib.pyplot as plt
import numpy as np
import pyart
from pyart.testing import get_test_data
warnings.filterwarnings('ignore')
## 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
Read in the Data and Plot the Data
Read in a sample file from the Southern Great Plains (SGP) Site
Our data is formatted as an mdv
file, which is the abbreviation for Meteorological Data Volume, a data format developed by the National Center for Atmospheric Research in the early 1990s.
file = get_test_data('110635.mdv')
radar = pyart.io.read(file)
Downloading file '110635.mdv' from 'https://adc.arm.gov/pyart/example_data/110635.mdv' to '/home/runner/.cache/pyart-datasets'.
Plot a quick-look of reflectivity and velocity
We can start by taking a quick look at the reflectivity and velocity fields. Notice how the velocity field is rather messy, indicated by the speckles and high/low values directly next to each other
fig = plt.figure(figsize=[8, 10])
ax = plt.subplot(211, projection=ccrs.PlateCarree())
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('reflectivity',
ax=ax,
sweep=0,
resolution='50m',
vmin=0,
vmax=60,
projection=ccrs.PlateCarree(),
cmap='pyart_HomeyerRainbow')
ax2 = plt.subplot(2,1,2,projection=ccrs.PlateCarree())
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('velocity',
ax=ax2,
sweep=0,
resolution='50m',
vmin=-17,
vmax=17,
projection=ccrs.PlateCarree(),
cmap='pyart_balance')
plt.show()
Dealiasing our Velocity
An Overview of Aliasing
The radial velocity measured by the radar is mesasured by detecing the phase shift between the transmitted pulse and the pulse recieved by the radar. However, using this methodology, it is only possible to detect phase shifts within \(\pm 2\pi\) due to the periodicity of the transmitted wave. Therefore, for example, a phase shift of \(3\pi\) would erroneously be detected as a phase shift of \(-\pi\) and give an incorrect value of velocity when retrieved by the radar. This phenomena is called aliasing. The maximium unambious velocity that can be detected by the radar before aliasing occurs is called the Nyquist velocity.
In the next example, you will see an example of aliasing occurring, where the values of +15 m/s abruptly transition into a region of -15 m/s, with -5 m/s in the middle of the region around 37 N, 97 W.
Calculate Velocity Texture
First, for dealiasing to work efficiently, we need to use a GateFilter. Notice that, this time, the data shown does not give us a nice gate_id. This is what raw data looks like, and we need to do some preprocessing on the data to remove noise and clutter. Thankfully, Py-ART gives us the capability to do this. As a simple filter in this example, we will first calculate the velocity texture using Py-ART’s calculate_velocity_texture
function. The velocity texture is the standard deviation of velocity surrounding a gate. This will be higher in the presence of noise.
Let’s investigate this function first…
pyart.retrieve.calculate_velocity_texture?
Determining the Right Parameters
You’ll notice that we need:
Our radar object
The name of our velocity field
The number of gates within our window to use to calculate the texture
The nyquist velocity
We can retrieve the nyquest velocity from our instrument parameters fortunately - using the following syntax!
nyquist_velocity = radar.instrument_parameters["nyquist_velocity"]["data"]
nyquist_velocity
array([16.52468, 16.52468, 16.52468, ..., 16.52468, 16.52468, 16.52468],
dtype=float32)
While the nyquist velocity is stored as an array, we see that these are all the same value…
np.unique(nyquist_velocity)
array([16.52468], dtype=float32)
Let’s save this single value to a float to use later…
nyquist_value = np.unique(nyquist_velocity)[0]
Calculate Velocity Texture and Filter our Data
Now that we have an ide?a of which parameters to pass in, let’s calculate velocity texture!
vel_texture = pyart.retrieve.calculate_velocity_texture(radar,
vel_field='velocity',
nyq=nyquist_value)
vel_texture
{'units': 'meters_per_second',
'standard_name': 'texture_of_radial_velocity_of_scatters_away_from_instrument',
'long_name': 'Doppler velocity texture',
'coordinates': 'elevation azimuth range',
'data': array([[5.8774082 , 5.73538228, 1.29254402, ..., 3.09248932, 3.11758189,
3.11758189],
[5.8774082 , 5.60215107, 1.29254402, ..., 7.30583692, 7.31494978,
7.3283932 ],
[6.98970645, 5.25008601, 1.93385618, ..., 7.32173353, 7.33318301,
7.33770852],
...,
[8.27919632, 7.05635664, 5.84890822, ..., 2.72774608, 2.76027457,
2.79050846],
[7.15030826, 6.57482555, 5.37450144, ..., 1.37623459, 1.40499757,
1.42364895],
[7.05635664, 6.57482555, 5.37450144, ..., 1.35595468, 1.35852404,
1.40499757]])}
The pyart.retrieve.calculate_velocity_texture
function results in a data dictionary, including the actual data, as well as metadata. We can add this to our radar
object, by using the radar.add_field
method, passing the name of our new field ("texture"
), the data dictionary (vel_texture
), and a clarification that we want to replace the existing velocity texture field if it already exists in our radar object (replace_existing=True
)
radar.add_field('texture', vel_texture, replace_existing=True)
Plot our Velocity Texture Field
Now that we have our velocity texture field added to our radar object, let’s plot it!
fig = plt.figure(figsize=[8, 8])
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('texture',
sweep=0,
resolution='50m',
vmin=0,
vmax=10,
projection=ccrs.PlateCarree(),
cmap='pyart_balance')
plt.show()
Determine a Suitable Velocity Texture Threshold
Plot a histogram of velocity texture to get a better idea of what values correspond to hydrometeors and what values of texture correspond to artifacts.
In the below example, a threshold of 3 would eliminate most of the peak corresponding to noise around 6 while preserving most of the values in the peak of ~0.5 corresponding to hydrometeors.
hist, bins = np.histogram(radar.fields['texture']['data'],
bins=150)
bins = (bins[1:]+bins[:-1])/2.0
plt.plot(bins,
hist,
label='Velocity Texture Frequency')
plt.axvline(3,
color='r',
label='Proposed Velocity Texture Threshold')
plt.xlabel('Velocity texture')
plt.ylabel('Count')
plt.legend()
<matplotlib.legend.Legend at 0x7fd6ae0210d0>
Setup a Gatefilter Object and Apply our Threshold
Now we can set up our GateFilter (pyart.filters.GateFilter
), which allows us to easily apply masks and filters to our radar object.
gatefilter = pyart.filters.GateFilter(radar)
gatefilter
<pyart.filters.gatefilter.GateFilter at 0x7fd6b6226650>
We discovered that a velocity texture threshold of only including values below 3 would be suitable for this dataset, we use the .exclude_above
method, specifying we want to exclude texture
values above 3.
gatefilter.exclude_above('texture', 3)
Plot our Filtered Data
Now that we have created a gatefilter, filtering our data using the velocity texture, let’s plot our data!
We need to pass our gatefilter to the plot_ppi_map
to apply it to our plot.
# Plot our Unfiltered Data
fig = plt.figure(figsize=[8, 10])
ax = plt.subplot(211, projection=ccrs.PlateCarree())
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('velocity',
title='Raw Radial Velocity (no filter)',
ax=ax,
sweep=0,
resolution='50m',
vmin=-17,
vmax=17,
projection=ccrs.PlateCarree(),
colorbar_label='Radial Velocity (m/s)',
cmap='pyart_balance')
ax2 = plt.subplot(2,1,2,projection=ccrs.PlateCarree())
# Plot our filtered data
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('velocity',
title='Radial Velocity with Velocity Texture Filter',
ax=ax2,
sweep=0,
resolution='50m',
vmin=-17,
vmax=17,
projection=ccrs.PlateCarree(),
colorbar_label='Radial Velocity (m/s)',
gatefilter=gatefilter,
cmap='pyart_balance')
plt.show()
Dealias the Velocity Using the Region-Based Method
At this point, we can use the dealias_region_based
to dealias the velocities and then add the new field to the radar!
velocity_dealiased = pyart.correct.dealias_region_based(radar,
vel_field='velocity',
nyquist_vel=nyquist_value,
centered=True,
gatefilter=gatefilter)
# Add our data dictionary to the radar object
radar.add_field('corrected_velocity', velocity_dealiased, replace_existing=True)
Plot our Cleaned, Dealiased Velocities
Plot the new velocities, which now look much more realistic.
fig = plt.figure(figsize=[8, 8])
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('corrected_velocity',
sweep=0,
resolution='50m',
vmin=-30,
vmax=30,
projection=ccrs.PlateCarree(),
colorbar_label='Radial Velocity (m/s)',
cmap='pyart_balance',
gatefilter=gatefilter)
plt.show()
Compare our Raw Velocity Field to our Dealiased, Cleaned Velocity Field
As a last comparison, let’s compare our raw, uncorrected velocities with our cleaned velocities, after applying the velocity texture threshold and dealiasing algorithm
# Plot our Unfiltered Data
fig = plt.figure(figsize=[8, 10])
ax = plt.subplot(211, projection=ccrs.PlateCarree())
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('velocity',
title='Raw Radial Velocity (no filter)',
ax=ax,
sweep=0,
resolution='50m',
vmin=-30,
vmax=30,
projection=ccrs.PlateCarree(),
colorbar_label='Radial Velocity (m/s)',
cmap='pyart_balance')
ax2 = plt.subplot(2,1,2,projection=ccrs.PlateCarree())
# Plot our filtered, dealiased data
display = pyart.graph.RadarMapDisplay(radar)
display.plot_ppi_map('corrected_velocity',
title='Radial Velocity with Velocity Texture Filter and Dealiasing',
ax=ax2,
sweep=0,
resolution='50m',
vmin=-30,
vmax=30,
projection=ccrs.PlateCarree(),
gatefilter=gatefilter,
colorbar_label='Radial Velocity (m/s)',
cmap='pyart_balance')
plt.show()
Conclusions
Within this lesson, we walked through how to apply radial velocity corrections to a dataset, filtering based on the velocity texture and using a regional dealiasing algorithm.
What’s Next
In the next few notebooks, we walk through gridding radar data and advanced visualization methods!
Resources and References
Py-ART essentials links: