Utilities (pyrad.util
)¶
Functions to read and write data and configuration files.
Radar Utilities¶
get_data_along_rng (radar, field_name, …[, …]) |
Get data at particular (azimuths, elevations) |
get_data_along_azi (radar, field_name, …[, …]) |
Get data at particular (ranges, elevations) |
get_data_along_ele (radar, field_name, …[, …]) |
Get data at particular (ranges, azimuths) |
get_ROI (radar, fieldname, sector) |
filter out any data outside the region of interest defined by sector |
rainfall_accumulation (t_in_vec, val_in_vec) |
Computes the rainfall accumulation of a time series over a given period |
time_series_statistics (t_in_vec, val_in_vec) |
Computes statistics over a time-averaged series. |
find_contiguous_times (times[, step]) |
Given and array of ordered times, find those contiguous according to a maximum time step |
join_time_series (t1, val1, t2, val2[, dropnan]) |
joins time_series. |
get_range_bins_to_avg (rad1_rng, rad2_rng) |
Compares the resolution of two radars and determines if and which radar has to be averaged and the length of the averaging window |
find_ray_index (ele_vec, azi_vec, ele, azi[, …]) |
Find the ray index corresponding to a particular elevation and azimuth |
find_rng_index (rng_vec, rng[, rng_tol]) |
Find the range index corresponding to a particular range |
find_nearest_gate (radar, lat, lon[, latlon_tol]) |
Find the radar gate closest to a lat,lon point |
find_neighbour_gates (radar, azi, rng[, …]) |
Find the neighbouring gates within +-delta_azi and +-delta_rng |
find_colocated_indexes (radar1, radar2, …) |
Given the theoretical elevation, azimuth and range of the co-located gates of two radars and a given tolerance returns the indices of the gates for the current radars |
get_target_elevations (radar_in) |
Gets RHI target elevations |
get_fixed_rng_data (radar, field_names, fixed_rng) |
Creates a 2D-grid with (azi, ele) data at a fixed range |
time_avg_range (timeinfo, avg_starttime, …) |
finds the new start and end time of an averaging |
get_closest_solar_flux (hit_datetime_list, …) |
finds the solar flux measurement closest to the sun hit |
create_sun_hits_field (rad_el, rad_az, …) |
creates a sun hits field from the position and power of the sun hits |
create_sun_retrieval_field (par, field_name, …) |
creates a sun retrieval field from the retrieval parameters |
compute_quantiles (field[, quantiles]) |
computes quantiles |
compute_quantiles_from_hist (bin_centers, hist) |
computes quantiles from histograms |
compute_quantiles_sweep (field, ray_start, …) |
computes quantiles of a particular sweep |
compute_2d_hist (field1, field2, field_name1, …) |
computes a 2D histogram of the data |
compute_1d_stats (field1, field2) |
returns statistics of data |
compute_2d_stats (field1, field2, …[, …]) |
computes a 2D histogram and statistics of the data |
compute_histogram (field, field_name[, …]) |
computes histogram of the data |
compute_histogram_sweep (field, ray_start, …) |
computes histogram of the data in a particular sweep |
belongs_roi_indices (lat, lon, roi) |
Get the indices of points that belong to roi in a list of points |
compute_profile_stats (field, gate_altitude, …) |
Compute statistics of vertical profile |
compute_directional_stats (field[, avg_type, …]) |
Computes the mean or the median along one of the axis (ray or range) |
project_to_vertical (data_in, data_height, …) |
Projects radar data to a regular vertical grid |
quantiles_weighted (values[, weight_vector, …]) |
Given a set of values and weights, compute the weighted quantile(s) and average. |
ratio_bootstrapping (nominator, denominator) |
Computes a set of samples obtained as sum(nominator)/sum(denominator) where the nominator and the denominator are randomly sampled with replacement. |
-
pyrad.util.
belongs_roi_indices
(lat, lon, roi)[source]¶ Get the indices of points that belong to roi in a list of points
Parameters: - lat, lon : float arrays
latitudes and longitudes to check
- roi : dict
Dictionary describing the region of interest
Returns: - inds : array of ints
list of indices of points belonging to ROI
- is_roi : str
Whether the list of points is within the region of interest. Can be ‘All’, ‘None’, ‘Some’
-
pyrad.util.
compute_1d_stats
(field1, field2)[source]¶ returns statistics of data
Parameters: - field1, field2 : ndarray 1D
the two fields to compare
Returns: - stats : dict
a dictionary with statistics
-
pyrad.util.
compute_2d_hist
(field1, field2, field_name1, field_name2, step1=None, step2=None)[source]¶ computes a 2D histogram of the data
Parameters: - field1, field2 : ndarray 2D
the radar fields
- field_name1, field_name2 : str
field names
- step1, step2 : float
size of the bins
Returns: - H : float array 2D
The bi-dimensional histogram of samples x and y
- xedges, yedges : float array
the bin edges along each dimension
-
pyrad.util.
compute_2d_stats
(field1, field2, field_name1, field_name2, step1=None, step2=None)[source]¶ computes a 2D histogram and statistics of the data
Parameters: - field1, field2 : ndarray 2D
the two fields
- field_name1, field_nam2: str
the name of the fields
- step1, step2 : float
size of bin
Returns: - hist_2d : array
the histogram
- bin_edges1, bin_edges2 : float array
The bin edges
- stats : dict
a dictionary with statistics
-
pyrad.util.
compute_directional_stats
(field, avg_type='mean', nvalid_min=1, axis=0)[source]¶ Computes the mean or the median along one of the axis (ray or range)
Parameters: - field : ndarray
the radar field
- avg_type :str
the type of average: ‘mean’ or ‘median’
- nvalid_min : int
the minimum number of points to consider the stats valid. Default 1
- axis : int
the axis along which to compute (0=ray, 1=range)
Returns: - values : ndarray 1D
The resultant statistics
- nvalid : ndarray 1D
The number of valid points used in the computation
-
pyrad.util.
compute_histogram
(field, field_name, bin_edges=None, step=None, vmin=None, vmax=None)[source]¶ computes histogram of the data
Parameters: - field : ndarray 2D
the radar field
- field_name: str or none
name of the field
- bins_edges :ndarray 1D
the bin edges
- step : float
size of bin
- vmin, vmax : float
The minimum and maximum value of the histogram
Returns: - bin_edges : float array
interval of each bin
- values : float array
values at each bin
-
pyrad.util.
compute_histogram_sweep
(field, ray_start, ray_end, field_name, step=None, vmin=None, vmax=None)[source]¶ computes histogram of the data in a particular sweep
Parameters: - field : ndarray 2D
the radar field
- ray_start, ray_end : int
starting and ending ray indexes
- field_name: str
name of the field
- step : float
size of bin
- vmin, vmax : float
minimum and maximum values
Returns: - bin_edges : float array
interval of each bin
- values : float array
values at each bin
-
pyrad.util.
compute_profile_stats
(field, gate_altitude, h_vec, h_res, quantity='quantiles', quantiles=array([0.25, 0.5, 0.75]), nvalid_min=4, std_field=None, np_field=None, make_linear=False, include_nans=False)[source]¶ Compute statistics of vertical profile
Parameters: - field : ndarray
the radar field
- gate_altitude: ndarray
the altitude at each radar gate [m MSL]
- h_vec : 1D ndarray
height vector [m MSL]
- h_res : float
heigh resolution [m]
- quantity : str
The quantity to compute. Can be [‘quantiles’, ‘mode’, ‘regression_mean’, ‘mean’]. If ‘mean’, the min, max, and average is computed.
- quantiles : 1D ndarray
the quantiles to compute
- nvalid_min : int
the minimum number of points to consider the stats valid
- std_field : ndarray
the standard deviation of the regression at each range gate
- np_field : ndarray
the number of points used to compute the regression at each range gate
- make_linear : Boolean
If true the data is transformed into linear coordinates before taking the mean
- include_nans : Boolean
If true NaN will be considered as zeros
Returns: - vals : ndarray 2D
The resultant statistics
- val_valid : ndarray 1D
The number of points to compute the stats used at each height level
-
pyrad.util.
compute_quantiles
(field, quantiles=None)[source]¶ computes quantiles
Parameters: - field : ndarray 2D
the radar field
- ray_start, ray_end : int
starting and ending ray indexes
- quantiles: float array
list of quantiles to compute
Returns: - quantiles : float array
list of quantiles
- values : float array
values at each quantile
-
pyrad.util.
compute_quantiles_from_hist
(bin_centers, hist, quantiles=None)[source]¶ computes quantiles from histograms
Parameters: - bin_centers : ndarray 1D
the bins
- hist : ndarray 1D
the histogram
- quantiles: float array
list of quantiles to compute
Returns: - quantiles : float array
list of quantiles
- values : float array
values at each quantile
-
pyrad.util.
compute_quantiles_sweep
(field, ray_start, ray_end, quantiles=None)[source]¶ computes quantiles of a particular sweep
Parameters: - field : ndarray 2D
the radar field
- ray_start, ray_end : int
starting and ending ray indexes
- quantiles: float array
list of quantiles to compute
Returns: - quantiles : float array
list of quantiles
- values : float array
values at each quantile
-
pyrad.util.
create_sun_hits_field
(rad_el, rad_az, sun_el, sun_az, data, imgcfg)[source]¶ creates a sun hits field from the position and power of the sun hits
Parameters: - rad_el, rad_az, sun_el, sun_az : ndarray 1D
azimuth and elevation of the radar and the sun respectively in degree
- data : masked ndarray 1D
the sun hit data
- imgcfg: dict
a dictionary specifying the ranges and resolution of the field to create
Returns: - field : masked ndarray 2D
the sun hit field
-
pyrad.util.
create_sun_retrieval_field
(par, field_name, imgcfg, lant=0.0)[source]¶ creates a sun retrieval field from the retrieval parameters
Parameters: - par : ndarray 1D
the 5 retrieval parameters
- imgcfg: dict
a dictionary specifying the ranges and resolution of the field to create
Returns: - field : masked ndarray 2D
the sun retrieval field
-
pyrad.util.
find_colocated_indexes
(radar1, radar2, rad1_ele, rad1_azi, rad1_rng, rad2_ele, rad2_azi, rad2_rng, ele_tol=0.5, azi_tol=0.5, rng_tol=50.0)[source]¶ Given the theoretical elevation, azimuth and range of the co-located gates of two radars and a given tolerance returns the indices of the gates for the current radars
Parameters: - radar1, radar2 : radar objects
the two radar objects
- rad1_ele, rad1_azi, rad1_rng : array of floats
the radar coordinates of the radar1 gates
- rad2_ele, rad2_azi, rad2_rng : array of floats
the radar coordinates of the radar2 gates
- ele_tol, azi_tol : floats
azimuth and elevation angle tolerance [deg]
- rng_tol : float
range Tolerance [m]
Returns: - ind_ray_rad1, ind_rng_rad1, ind_ray_rad2, ind_rng_rad2 : array of ints
the ray and range indexes of each radar gate
-
pyrad.util.
find_contiguous_times
(times, step=600)[source]¶ Given and array of ordered times, find those contiguous according to a maximum time step
Parameters: - times : array of datetimes
The array of times
- step : float
The time step [s]
Returns: - start_times, end_times : array of date times
The start and end of each consecutive time period
-
pyrad.util.
find_nearest_gate
(radar, lat, lon, latlon_tol=0.0005)[source]¶ Find the radar gate closest to a lat,lon point
Parameters: - radar : radar object
the radar object
- lat, lon : float
The position of the point
- latlon_tol : float
The tolerance around this point
Returns: - ind_ray, ind_rng : int
The ray and range index
- azi, rng : float
the range and azimuth position of the gate
-
pyrad.util.
find_neighbour_gates
(radar, azi, rng, delta_azi=None, delta_rng=None)[source]¶ Find the neighbouring gates within +-delta_azi and +-delta_rng
Parameters: - radar : radar object
the radar object
- azi, rng : float
The azimuth [deg] and range [m] of the central gate
- delta_azi, delta_rng : float
The extend where to look for
Returns: - inds_ray_aux, ind_rng_aux : int
The indices (ray, rng) of the neighbouring gates
-
pyrad.util.
find_ray_index
(ele_vec, azi_vec, ele, azi, ele_tol=0.0, azi_tol=0.0, nearest='azi')[source]¶ Find the ray index corresponding to a particular elevation and azimuth
Parameters: - ele_vec, azi_vec : float arrays
The elevation and azimuth data arrays where to look for
- ele, azi : floats
The elevation and azimuth to search
- ele_tol, azi_tol : floats
Tolerances [deg]
- nearest : str
criteria to define wich ray to keep if multiple rays are within tolerance. azi: nearest azimuth, ele: nearest elevation
Returns: - ind_ray : int
The ray index
-
pyrad.util.
find_rng_index
(rng_vec, rng, rng_tol=0.0)[source]¶ Find the range index corresponding to a particular range
Parameters: - rng_vec : float array
The range data array where to look for
- rng : float
The range to search
- rng_tol : float
Tolerance [m]
Returns: - ind_rng : int
The range index
-
pyrad.util.
get_ROI
(radar, fieldname, sector)[source]¶ filter out any data outside the region of interest defined by sector
Parameters: - radar : radar object
the radar object where the data is
- fieldname : str
name of the field to filter
- sector : dict
a dictionary defining the region of interest
Returns: - roi_flag : ndarray
a field array with ones in gates that are in the Region of Interest
-
pyrad.util.
get_closest_solar_flux
(hit_datetime_list, flux_datetime_list, flux_value_list)[source]¶ finds the solar flux measurement closest to the sun hit
Parameters: - hit_datetime_list : datetime array
the date and time of the sun hit
- flux_datetime_list : datetime array
the date and time of the solar flux measurement
- flux_value_list: ndarray 1D
the solar flux values
Returns: - flux_datetime_closest_list : datetime array
the date and time of the solar flux measurement closest to sun hit
- flux_value_closest_list : ndarray 1D
the solar flux values closest to the sun hit time
-
pyrad.util.
get_data_along_azi
(radar, field_name, fix_ranges, fix_elevations, rng_tol=50.0, ang_tol=1.0, azi_start=None, azi_stop=None)[source]¶ Get data at particular (ranges, elevations)
Parameters: - radar : radar object
the radar object where the data is
- field_name : str
name of the field to filter
- fix_ranges, fix_elevations: list of floats
List of ranges [m], elevations [deg] couples
- rng_tol : float
Tolerance between the nominal range and the radar range [m]
- ang_tol : float
Tolerance between the nominal angle and the radar angle [deg]
- azi_start, azi_stop: float
Start and stop azimuth angle of the data [deg]
Returns: - xvals : list of float arrays
The ranges of each rng, ele pair
- yvals : list of float arrays
The values
- valid_rng, valid_ele : float arrays
The rng, ele pairs
-
pyrad.util.
get_data_along_ele
(radar, field_name, fix_ranges, fix_azimuths, rng_tol=50.0, ang_tol=1.0, ele_min=None, ele_max=None)[source]¶ Get data at particular (ranges, azimuths)
Parameters: - radar : radar object
the radar object where the data is
- field_name : str
name of the field to filter
- fix_ranges, fix_azimuths: list of floats
List of ranges [m], azimuths [deg] couples
- rng_tol : float
Tolerance between the nominal range and the radar range [m]
- ang_tol : float
Tolerance between the nominal angle and the radar angle [deg]
- ele_min, ele_max: float
Min and max elevation angle [deg]
Returns: - xvals : list of float arrays
The ranges of each rng, ele pair
- yvals : list of float arrays
The values
- valid_rng, valid_ele : float arrays
The rng, ele pairs
-
pyrad.util.
get_data_along_rng
(radar, field_name, fix_elevations, fix_azimuths, ang_tol=1.0, rmin=None, rmax=None)[source]¶ Get data at particular (azimuths, elevations)
Parameters: - radar : radar object
the radar object where the data is
- field_name : str
name of the field to filter
- fix_elevations, fix_azimuths: list of floats
List of elevations, azimuths couples [deg]
- ang_tol : float
Tolerance between the nominal angle and the radar angle [deg]
- rmin, rmax: float
Min and Max range of the obtained data [m]
Returns: - xvals : list of float arrays
The ranges of each azi, ele pair
- yvals : list of float arrays
The values
- valid_azi, valid_ele : float arrays
The azi, ele pairs
-
pyrad.util.
get_fixed_rng_data
(radar, field_names, fixed_rng, rng_tol=50.0, ele_min=None, ele_max=None, azi_min=None, azi_max=None)[source]¶ Creates a 2D-grid with (azi, ele) data at a fixed range
Parameters: - radar : radar object
The radar object containing the data
- field_name : str
The field name
- fixed_rng : float
The fixed range [m]
- rng_tol : float
The tolerance between the nominal range and the actual radar range [m]
- 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
-
pyrad.util.
get_range_bins_to_avg
(rad1_rng, rad2_rng)[source]¶ Compares the resolution of two radars and determines if and which radar has to be averaged and the length of the averaging window
Parameters: - rad1_rng : array
the range of radar 1
- rad2_rng : datetime
the range of radar 2
Returns: - avg_rad1, avg_rad2 : Boolean
Booleans specifying if the radar data has to be average in range
- avg_rad_lim : array with two elements
the limits to the average (centered on each range gate)
-
pyrad.util.
get_target_elevations
(radar_in)[source]¶ Gets RHI target elevations
Parameters: - radar_in : Radar object
current radar object
Returns: - target_elevations : 1D-array
Azimuth angles
- el_tol : float
azimuth tolerance
-
pyrad.util.
join_time_series
(t1, val1, t2, val2, dropnan=False)[source]¶ joins time_series. Only of package pandas is available otherwise returns None.
Parameters: - t1 : datetime array
time of first series
- val1 : float array
value of first series
- t2 : datetime array
time of second series
- val2 : float array
value of second series
- dropnan : boolean
if True remove NaN from the time series
Returns: - t_out_vec : datetime array
the resultant date time after joining the series
- val1_out_vec : float array
value of first series
- val2_out_vec : float array
value of second series
-
pyrad.util.
project_to_vertical
(data_in, data_height, grid_height, interp_kind='none', fill_value=-9999.0)[source]¶ Projects radar data to a regular vertical grid
Parameters: - data_in : ndarray 1D
the radar data to project
- data_height : ndarray 1D
the height of each radar point
- grid_height : ndarray 1D
the regular vertical grid to project to
- interp_kind : str
The type of interpolation to use: ‘none’ or ‘nearest’
- fill_value : float
The fill value used for interpolation
Returns: - data_out : ndarray 1D
The projected data
-
pyrad.util.
quantiles_weighted
(values, weight_vector=None, quantiles=array([0.5]), weight_threshold=None, data_is_log=False, nvalid_min=3)[source]¶ Given a set of values and weights, compute the weighted quantile(s) and average.
Parameters: - values : array of floats
Array containing the values. Can be 2-dimensional
- weight_vector : array of floats or None
array containing the weights to apply. If None it will be an array of ones (uniform weight). If values is a 2D array it will be repeated for the second dimension
- quantiles : array of floats
The quantiles to be computed
- weight_threshold : float or None
If weight_threshold is set quantiles will be computed only if the total weight (sum of the weights of valid data) exceeds this threshold
- data_is_log : Bool
If true the values will be considered to be in logarithmic scale and transformed into linear scale before computing the quantiles and average
- nvalid_min : int
Minimum number of valid points to consider the computation valid
Returns: - avg : float
the weighted average
- quants : array of floats
an array containing the weighted quantiles in the same order as the quantiles vector
- nvalid : int
Number of valid points in the computation of the statistics
-
pyrad.util.
rainfall_accumulation
(t_in_vec, val_in_vec, cum_time=3600.0, base_time=0.0, dropnan=False)[source]¶ Computes the rainfall accumulation of a time series over a given period
Parameters: - t_in_vec : datetime array
the input date and time array
- val_in_vec : float array
the input values array [mm/h]
- cum_time : int
accumulation time [s]
- base_time : int
base time [s]
- dropnan : boolean
if True remove NaN from the time series
Returns: - t_out_vec : datetime array
the output date and time array
- val_out_vec : float array
the output values array
- np_vec : int array
the number of samples at each period
-
pyrad.util.
ratio_bootstrapping
(nominator, denominator, nsamples=1000)[source]¶ Computes a set of samples obtained as sum(nominator)/sum(denominator) where the nominator and the denominator are randomly sampled with replacement.
Parameters: - nominator, denominator : 1D array
The data points in the nominator and the denominator. Nominator and denominator are not independent, i.e. data point i in the nominator is linked to data point i in the denominator
- nsamples : int
Number of iteration, i.e. number of samples desired
Returns: - samples : 1D array
the resultant samples
-
pyrad.util.
time_avg_range
(timeinfo, avg_starttime, avg_endtime, period)[source]¶ finds the new start and end time of an averaging
Parameters: - timeinfo : datetime
the current volume time
- avg_starttime : datetime
the current average start time
- avg_endtime: datetime
the current average end time
- period: float
the averaging period
Returns: - new_starttime : datetime
the new average start time
- new_endtime : datetime
the new average end time
-
pyrad.util.
time_series_statistics
(t_in_vec, val_in_vec, avg_time=3600, base_time=1800, method='mean', dropnan=False)[source]¶ Computes statistics over a time-averaged series. Only of package pandas is available otherwise returns None
Parameters: - t_in_vec : datetime array
the input date and time array
- val_in_vec : float array
the input values array
- avg_time : int
averaging time [s]
- base_time : int
base time [s]
- method : str
statistical method
- dropnan : boolean
if True remove NaN from the time series
Returns: - t_out_vec : datetime array
the output date and time array
- val_out_vec : float array
the output values array