Utilities (pyart.util)

Miscellaneous utility functions.

The location and names of these functions within Py-ART may change between versions without depreciation, use with caution.

Direction statistics

angular_mean(angles) Compute the mean of a distribution of angles in radians.
angular_std(angles) Compute the standard deviation of a distribution of angles in radians.
angular_mean_deg(angles) Compute the mean of a distribution of angles in degrees.
angular_std_deg(angles) Compute the standard deviation of a distribution of angles in degrees.
interval_mean(dist, interval_min, interval_max) Compute the mean of a distribution within an interval.
interval_std(dist, interval_min, interval_max) Compute the standard deviation of a distribution within an interval.
mean_of_two_angles(angles1, angles2) Compute the element by element mean of two sets of angles.
mean_of_two_angles_deg(angle1, angle2) Compute the element by element mean of two sets of angles in degrees.

Miscellaneous functions

cross_section_ppi(radar, target_azimuths[, …]) Extract cross sections from a PPI volume along one or more azimuth angles.
cross_section_rhi(radar, target_elevations) Extract cross sections from an RHI volume along one or more elevation angles.
colocated_gates(radar1, radar2[, h_tol, …]) Flags radar gates of radar1 colocated with radar2
intersection(radar1, radar2[, h_tol, …]) Flags region of radar1 that is intersecting with radar2 and complies with criteria regarding visibility, altitude, range, elevation angle and azimuth angle
datetime_from_radar(radar[, epoch]) Return a datetime for the first ray in a Radar.
datetimes_from_radar(radar[, epoch]) Return an array of datetimes for the rays in a Radar.
datetime_from_dataset(dataset[, epoch]) Return a datetime for the first time in a netCDF Dataset.
datetimes_from_dataset(dataset[, epoch]) Return an array of datetimes for the times in a netCDF Dataset.
datetime_from_grid(grid[, epoch]) Return a datetime for the volume start in a Grid.
estimate_noise_hs74(spectrum[, navg, nnoise_min]) Estimate noise parameters of a Doppler spectrum.
estimate_noise_ivic13(pwr_w_ray[, pct, …]) Estimate noise parameters of a ray
get_ivic_pct(npulses, pct, z[, prob_thr]) Get the point clutter threshold (PCT)
get_ivic_flat_reg_var_max(npulses, y, x, …) Get the threshold for maximum local variance of noise [dB]
get_ivic_snr_thr(npulses, snr_thr[, pfa_thr]) Get the threshold for steps 3 and 6 of ivic
ivic_pct_table(npulses) Get the point clutter threshold (PCT) of Ivic from a look up table.
ivic_snr_thr_table(npulses) Get the threshold for steps 3 and 6 of ivic from a look up table The thresholds are computed for between 3 and 200 pulses.
ivic_flat_reg_var_max_table(npulses, …) Get maximum variance of noise of Ivic from a look up table.
is_vpt(radar[, offset]) Determine if a Radar appears to be a vertical pointing scan.
to_vpt(radar[, single_scan]) Convert an existing Radar object to represent a vertical pointing scan.
join_radar(radar1, radar2) Combine two radar instances into one.
join_spectra(spectra1, spectra2) Combine two spectra instances into one.
cut_radar(radar, field_names[, rng_min, …]) Cuts the radar object into new dimensions
cut_radar_spectra(radar, field_names[, …]) Cuts the radar spectra object into new dimensions
radar_from_spectra(psr) obtain a Radar object from a RadarSpectra object
interpol_spectra(psr[, kind, fill_value]) Interpolates the spectra so that it has a uniform grid
simulated_vel_from_profile(radar, profile[, …]) Create simulated radial velocities from a profile of horizontal winds.
texture_along_ray(radar, var[, wind_size]) Compute field texture along ray using a user specified window size.
texture(radar, var) Determine a texture field using an 11pt stdev texarray=texture(pyradarobj, field).
rolling_window(a, window) create a rolling window object for application of functions eg: result=np.ma.std(array, 11), 1)
angular_texture_2d(image, N, interval) Compute the angular texture of an image.
pyart.util.angular_mean(angles)[source]

Compute the mean of a distribution of angles in radians.

Parameters:
angles : array like

Distribution of angles in radians.

Returns:
mean : float

The mean angle of the distribution in radians.

pyart.util.angular_mean_deg(angles)[source]

Compute the mean of a distribution of angles in degrees.

Parameters:
angles : array like

Distribution of angles in degrees.

Returns:
mean : float

The mean angle of the distribution in degrees.

pyart.util.angular_std(angles)[source]

Compute the standard deviation of a distribution of angles in radians.

Parameters:
angles : array like

Distribution of angles in radians.

Returns:
std : float

Standard deviation of the distribution.

pyart.util.angular_std_deg(angles)[source]

Compute the standard deviation of a distribution of angles in degrees.

Parameters:
angles : array like

Distribution of angles in degrees.

Returns:
std : float

Standard deviation of the distribution.

pyart.util.angular_texture_2d(image, N, interval)[source]

Compute the angular texture of an image. Uses convolutions in order to speed up texture calculation by a factor of ~50 compared to using ndimage.generic_filter.

Parameters:
image : 2D array of floats

The array containing the velocities in which to calculate texture from.

N : int

This is the window size for calculating texture. The texture will be calculated from an N by N window centered around the gate.

interval : float

The absolute value of the maximum velocity. In conversion to radial coordinates, pi will be defined to be interval and -pi will be -interval. It is recommended that interval be set to the Nyquist velocity.

Returns:
std_dev : float array

Texture of the radial velocity field.

pyart.util.colocated_gates(radar1, radar2, h_tol=0.0, latlon_tol=0.0, coloc_gates_field=None)[source]

Flags radar gates of radar1 colocated with radar2

Parameters:
radar1 : Radar

radar object that is going to be flagged

radar2 : Radar

radar object

h_tol : float

tolerance in altitude [m]

latlon_tol : float

tolerance in latitude/longitude [deg]

coloc_gates_field : string

Name of the field to retrieve the data

Returns:
coloc_dict : dict

a dictionary containing the colocated positions of radar 1 (ele, azi, rng) and radar 2

coloc_rad1 :

field with the colocated gates of radar1 flagged, i.e: 1: not colocated gates 2: colocated (0 is reserved)

pyart.util.cross_section_ppi(radar, target_azimuths, az_tol=None)[source]

Extract cross sections from a PPI volume along one or more azimuth angles.

Parameters:
radar : Radar

Radar volume containing PPI sweeps from which azimuthal cross sections will be extracted.

target_azimuth : list

Azimuthal angles in degrees where cross sections will be taken.

az_tol : float, optional

Azimuth angle tolerance in degrees. If none the nearest angle is used. If valid only angles within the tolerance distance are considered.

Returns:
radar_rhi : Radar

Radar volume containing RHI sweeps which contain azimuthal cross sections from the original PPI volume.

pyart.util.cross_section_rhi(radar, target_elevations, el_tol=None)[source]

Extract cross sections from an RHI volume along one or more elevation angles.

Parameters:
radar : Radar

Radar volume containing RHI sweeps from which azimuthal cross sections will be extracted.

target_elevations : list

Elevation angles in degrees where cross sections will be taken.

el_tol : float, optional

Elevation angle tolerance in degrees. If none the nearest angle is used. If valid only angles within the tolerance distance are considered.

Returns:
radar_ppi : Radar

Radar volume containing PPI sweeps which contain azimuthal cross sections from the original RHI volume.

pyart.util.cut_radar(radar, field_names, rng_min=None, rng_max=None, ele_min=None, ele_max=None, azi_min=None, azi_max=None)[source]

Cuts the radar object into new dimensions

Parameters:
radar : radar object

The radar object containing the data

field_names : str or None

The fields to keep in the new radar

rng_min, rng_max : float

The range limits [m]. If None the entire coverage of the radar is going to be used

ele_min, ele_max, azi_min, azi_max : float or None

The limits of the grid [deg]. If None the limits will be the limits of the radar volume

Returns:
radar : radar object

The radar object containing only the desired data

pyart.util.cut_radar_spectra(radar, field_names, rng_min=None, rng_max=None, ele_min=None, ele_max=None, azi_min=None, azi_max=None)[source]

Cuts the radar spectra object into new dimensions

Parameters:
radar : radar object

The radar object containing the data

field_names : str or None

The fields to keep in the new radar

rng_min, rng_max : float

The range limits [m]. If None the entire coverage of the radar is going to be used

ele_min, ele_max, azi_min, azi_max : float or None

The limits of the grid [deg]. If None the limits will be the limits of the radar volume

Returns:
radar : radar object

The radar object containing only the desired data

pyart.util.datetime_from_dataset(dataset, epoch=False, **kwargs)[source]

Return a datetime for the first time in a netCDF Dataset.

pyart.util.datetime_from_grid(grid, epoch=False, **kwargs)[source]

Return a datetime for the volume start in a Grid.

pyart.util.datetime_from_radar(radar, epoch=False, **kwargs)[source]

Return a datetime for the first ray in a Radar.

pyart.util.datetimes_from_dataset(dataset, epoch=False, **kwargs)[source]

Return an array of datetimes for the times in a netCDF Dataset.

pyart.util.datetimes_from_radar(radar, epoch=False, **kwargs)[source]

Return an array of datetimes for the rays in a Radar.

pyart.util.estimate_noise_hs74(spectrum, navg=1, nnoise_min=1)[source]

Estimate noise parameters of a Doppler spectrum.

Use the method of estimating the noise level in Doppler spectra outlined by Hildebrand and Sehkon, 1974.

Parameters:
spectrum : array like

Doppler spectrum in linear units.

navg : int, optional

The number of spectral bins over which a moving average has been taken. Corresponds to the p variable from equation 9 of the article. The default value of 1 is appropiate when no moving average has been applied to the spectrum.

nnoise_min : int

Minimum number of noise samples to consider the estimation valid

Returns:
mean : float-like

Mean of points in the spectrum identified as noise.

threshold : float-like

Threshold separating noise from signal. The point in the spectrum with this value or below should be considered as noise, above this value signal. It is possible that all points in the spectrum are identified as noise. If a peak is required for moment calculation then the point with this value should be considered as signal.

var : float-like

Variance of the points in the spectrum identified as noise.

nnoise : int

Number of noise points in the spectrum.

References

P. H. Hildebrand and R. S. Sekhon, Objective Determination of the Noise Level in Doppler Spectra. Journal of Applied Meteorology, 1974, 13, 808-811.

pyart.util.estimate_noise_ivic13(pwr_w_ray, pct=2.0, delay=1, flat_reg_wlen=96, flat_reg_var_max=1.5, snr_thr=1.6, npulses=30, ngates_min=800, iterations=10, get_noise_pos=False)[source]

Estimate noise parameters of a ray

Use the method of estimating the noise level in a ray outlined by Ivic, 2013.

Parameters:
pwr_w_ray : array like

Doppler spectrum in linear units.

pct : float

Point Clutter Threshold

delay : int

distance of the gate to which compare for point target threshold

flat_reg_wlen : int

Minimum number of gates that should contain a valid region. Default gives a size of 8 km with 83.3 m resolution

flat_reg_var_max : float

Maximum local variance of powers in decibels to consider the region as flat.

snr_thr : float

Threshold applied in steps 3 and 6

npulses : int

Number of pulses used to compute the power of the array

ngates_min: int

minimum number of gates with noise to consider the retrieval valid

iterations: int

number of iterations in step 7

get_noise_pos : bool

If True the indices of the gates identified as noise will be returned

Returns:
mean : float-like

Mean of the gates in the ray identified as noise.

var : float-like

Variance of the gates in the ray identified as noise.

nnoise : int

Number of noise gates in the ray.

inds_ray : 1D-array

The indices of the gates containing noise

References

I.R. Ivic, C. Curtis and S.M. Torres, Radial-Based Noise Power Estimation for Weather Radars. Journal of Atmospheric and Oceanic Technology, 2013, 30, 2737-2753.

pyart.util.get_ivic_flat_reg_var_max(npulses, y, x, flat_reg_wlen, ndBm=30.0, accuracy=0.001)[source]

Get the threshold for maximum local variance of noise [dB]

Parameters:
npulses : 1D array

array with the number of pulses in a ray

y : 1D array

array of y values

x : 1D array

array of x values. These are the possible values of the threshold

flat_reg_wlen : int

the lenght of the flat region window in bins

ndBm : float

the mean noise power in dBm

accuracy : float

The desired accuracy of the threshold (difference with target probability)

Returns:
var_thr : 1D array

The thresholds corresponding to each number of pulses

pyart.util.get_ivic_pct(npulses, pct, z, prob_thr=0.0001)[source]

Get the point clutter threshold (PCT)

Parameters:
npulses : 1D array

array with the number of pulses in a ray

pct : 1D array

array with possible values of the PCT

z : 1D array

array of z values

prob_thr : float

Desired probability

Returns:
pct_out : 1D array

The PCT threshold corresponding to each number of pulses

pyart.util.get_ivic_snr_thr(npulses, snr_thr, pfa_thr=0.001)[source]

Get the threshold for steps 3 and 6 of ivic

Parameters:
npulses : 1D array

array with the number of pulses in a ray

snr_thr : 1D array

array with possible values of the snr threshold

pfa_thr : float

the desired probability of false alarm

Returns:
snr_thr_out : 1D array

The snr threshold corresponding to each number of pulses according to the desired probability of false alarm

pyart.util.interpol_spectra(psr, kind='linear', fill_value=0.0)[source]

Interpolates the spectra so that it has a uniform grid

Parameters:
psr : RadarSpectra object

The original spectra

Returns:
psr_interp : RadarSpectra object

The interpolated spectra

pyart.util.intersection(radar1, radar2, h_tol=0.0, latlon_tol=0.0, vol_d_tol=None, vismin=None, hmin=None, hmax=None, rmin=None, rmax=None, elmin=None, elmax=None, azmin=None, azmax=None, visib_field=None, intersec_field=None)[source]

Flags region of radar1 that is intersecting with radar2 and complies with criteria regarding visibility, altitude, range, elevation angle and azimuth angle

Parameters:
radar1 : Radar

radar object that is going to be flagged

radar2 : Radar

radar object checked for intersecting region

h_tol : float

tolerance in altitude [m]

latlon_tol : float

latitude and longitude tolerance [decimal deg]

vol_d_tol : float

pulse volume diameter tolerance [m]

vismin : float

minimum visibility [percentage]

hmin, hmax : floats

min and max altitude [m MSL]

rmin, rmax : floats

min and max range from radar [m]

elmin, elmax : floats

min and max elevation angle [deg]

azmin, azmax : floats

min and max azimuth angle [deg]

Returns:
intersec_rad1_dict : dict

the field with the gates of radar1 in the same region as radar2 flagged, i.e.: 1 not intersecting, 2 intersecting, 0 is reserved

pyart.util.interval_mean(dist, interval_min, interval_max)[source]

Compute the mean of a distribution within an interval.

Return the average of the array elements which are interpreted as being taken from a circular interval with endpoints given by interval_min and interval_max.

Parameters:
dist : array like

Distribution of values within an interval.

interval_min, interval_max : float

The endpoints of the interval.

Returns:
mean : float

The mean value of the distribution.

pyart.util.interval_std(dist, interval_min, interval_max)[source]

Compute the standard deviation of a distribution within an interval.

Return the standard deviation of the array elements which are interpreted as being taken from a circular interval with endpoints given by interval_min and interval_max.

Parameters:
dist : array_like

Distribution of values within an interval.

interval_min, interval_max : float

The endpoints of the interval.

Returns:
std : float

The standard deviation of the distribution.

pyart.util.is_vpt(radar, offset=0.5)[source]

Determine if a Radar appears to be a vertical pointing scan.

This function only verifies that the object is a vertical pointing scan, use the to_vpt() function to convert the radar to a vpt scan if this function returns True.

Parameters:
radar : Radar

Radar object to determine if.

offset : float, optional

Maximum offset of the elevation from 90 degrees to still consider to be vertically pointing.

Returns:
flag : bool

True if the radar appear to be verticle pointing, False if not.

pyart.util.ivic_flat_reg_var_max_table(npulses, flat_reg_wlen)[source]

Get maximum variance of noise of Ivic from a look up table.

Parameters:
npulses : 1D array

array with the number of pulses in a ray

flat_reg_wlen : float

The flat region window length in bins

Returns:
flat_reg_var_max : 1D array

The maximum variance threshold corresponding to each number of pulses

pyart.util.ivic_pct_table(npulses)[source]

Get the point clutter threshold (PCT) of Ivic from a look up table. The thresholds are computed for between 3 and 78 pulses. If there number of pulses is beyond this range it throws an error

Parameters:
npulses : 1D array

array with the number of pulses in a ray

Returns:
pct_table : 1D array

The PCT corresponding to each number of pulses

pyart.util.ivic_snr_thr_table(npulses)[source]

Get the threshold for steps 3 and 6 of ivic from a look up table The thresholds are computed for between 3 and 200 pulses. If there number of pulses is beyond this range it throws an error

Parameters:
npulses : 1D array

array with the number of pulses in a ray

Returns:
snr_thr_out : 1D array

The snr threshold corresponding to each number of pulses according to the desired probability of false alarm

pyart.util.join_radar(radar1, radar2)[source]

Combine two radar instances into one.

Parameters:
radar1 : Radar

Radar object.

radar2 : Radar

Radar object.

pyart.util.join_spectra(spectra1, spectra2)[source]

Combine two spectra instances into one.

Parameters:
spectra1 : spectra

spectra object.

spectra2 : spectra

spectra object.

pyart.util.mean_of_two_angles(angles1, angles2)[source]

Compute the element by element mean of two sets of angles.

Parameters:
angles1 : array

First set of angles in radians.

angles2 : array

Second set of angles in radians.

Returns:
mean : array

Elements by element angular mean of the two sets of angles in radians.

pyart.util.mean_of_two_angles_deg(angle1, angle2)[source]

Compute the element by element mean of two sets of angles in degrees.

Parameters:
angle1 : array

First set of angles in degrees.

angle2 : array

Second set of angles in degrees.

Returns:
mean : array

Elements by element angular mean of the two sets of angles in degrees.

pyart.util.radar_from_spectra(psr)[source]

obtain a Radar object from a RadarSpectra object

Parameters:
psr : RadarSpectra object

The reference object

Returns:
radar : radar object

The new radar object

pyart.util.rolling_window(a, window)[source]

create a rolling window object for application of functions eg: result=np.ma.std(array, 11), 1)

pyart.util.simulated_vel_from_profile(radar, profile, interp_kind='linear', sim_vel_field=None)[source]

Create simulated radial velocities from a profile of horizontal winds.

Parameters:
radar : Radar

Radar instance which provides the scanning parameters for the simulated radial velocities.

profile : HorizontalWindProfile

Profile of horizontal winds.

interp_kind : str, optional

Specifies the kind of interpolation used to determine the winds at a given height. Must be one of ‘linear’, ‘nearest’, ‘zero’, ‘slinear’, ‘quadratic’, or ‘cubic’. The the documentation for the SciPy scipy.interpolate.interp1d function for descriptions.

sim_vel_field : str, optional

Name to use for the simulated velocity field metadata. None will use the default field name from the Py-ART configuration file.

Returns:
sim_vel : dict

Dictionary containing a radar field of simulated radial velocities.

pyart.util.texture(radar, var)[source]

Determine a texture field using an 11pt stdev texarray=texture(pyradarobj, field).

pyart.util.texture_along_ray(radar, var, wind_size=7)[source]

Compute field texture along ray using a user specified window size.

Parameters:
radar : radar object

The radar object where the field is.

var : str

Name of the field which texture has to be computed.

wind_size : int, optional

Optional. Size of the rolling window used.

Returns:
tex : radar field

The texture of the specified field.

pyart.util.to_vpt(radar, single_scan=True)[source]

Convert an existing Radar object to represent a vertical pointing scan.

This function does not verify that the Radar object contains a vertical pointing scan. To perform such a check use is_vpt().

Parameters:
radar : Radar

Mislabeled vertical pointing scan Radar object to convert to be properly labeled. This object is converted in place, no copy of the existing data is made.

single_scan : bool, optional

True to convert the volume to a single scan, any azimuth angle data is lost. False will convert the scan to contain the same number of scans as rays, azimuth angles are retained.