Radar Corrections (pyart.correct)

Correct radar fields.

Velocity unfolding

dealias_fourdd(radar[, last_radar, …]) Dealias Doppler velocities using the 4DD algorithm.
dealias_unwrap_phase(radar[, unwrap_unit, …]) Dealias Doppler velocities using multi-dimensional phase unwrapping.
dealias_region_based(radar[, ref_vel_field, …]) Dealias Doppler velocities using a region based algorithm.

Other corrections

calculate_attenuation(radar, z_offset[, …]) Calculate the attenuation from a polarimetric radar using Z-PHI method.
calculate_attenuation_zphi(radar[, doc, …]) Calculate the attenuation and the differential attenuation from a polarimetric radar using Z-PHI method..
calculate_attenuation_philinear(radar[, …]) Calculate the attenuation and the differential attenuation from a polarimetric radar using linear dependece with PhiDP.
phase_proc_lp(radar, offset[, debug, …]) Phase process using a LP method [1].
det_sys_phase_ray(radar[, ind_rmin, …]) Public method Alternative determination of the system phase.
correct_sys_phase(radar[, ind_rmin, …]) correction of the system offset.
smooth_phidp_single_window(radar[, …]) correction of the system offset and smoothing using one window
smooth_phidp_double_window(radar[, …]) correction of the system offset and smoothing using two window
despeckle_field(radar, field[, label_dict, …]) Despeckle a radar volume by identifying small objects in each scan and masking them out.
correct_noise_rhohv(radar[, urhohv_field, …]) Corrects RhoHV for noise according to eq.
correct_bias(radar[, bias, field_name]) Corrects a radar data bias.
correct_visibility(radar[, vis_field, …]) Corrects the reflectivity according to visibility.
est_rhohv_rain(radar[, ind_rmin, ind_rmax, …]) Estimates the quantiles of RhoHV in rain for each sweep
est_zdr_precip(radar[, ind_rmin, ind_rmax, …]) Filters out all undesired data to be able to estimate ZDR bias, either in moderate rain or from vertically pointing scans
est_zdr_snow(radar[, ind_rmin, ind_rmax, …]) Filters out all undesired data to be able to estimate ZDR bias in snow
selfconsistency_bias(radar, zdr_kdpzh_dict) Estimates reflectivity bias at each ray using the self-consistency algorithm by Gourley
selfconsistency_bias2(radar, zdr_kdpzh_dict) Estimates reflectivity bias at each ray using the self-consistency algorithm by Gourley
selfconsistency_kdp_phidp(radar, zdr_kdpzh_dict) Estimates KDP and PhiDP in rain from Zh and ZDR using a selfconsistency relation between ZDR, Zh and KDP.
get_sun_hits(radar[, delev_max, dazim_max, …]) get data from suspected sun hits.
get_sun_hits_ivic(radar[, delev_max, …]) get data from suspected sun hits.
sun_retrieval(az_rad, az_sun, el_rad, …[, …]) Estimates sun parameters from sun hits
phase_proc_lp_gf(radar[, gatefilter, debug, …]) Phase process using a LP method [1] using Py-ART’s Gatefilter.

Helper functions

find_objects(radar, field, threshold[, …]) Find objects (i.e., contiguous gates) in one or more sweeps that match thresholds.
get_mask_fzl(radar[, fzl, doc, min_temp, …]) Constructs a mask to mask data placed thickness m below data at min_temp and beyond.
sun_power(solar_flux, pulse_width, wavelen, …) computes the theoretical sun power detected at the antenna [dBm] as it would be without atmospheric attenuation (sun power at top of the atmosphere) for a given solar flux and radar characteristics
ptoa_to_sf(ptoa, pulse_width, wavelen, …) Converts the sun power at the top of the atmosphere (in dBm) into solar flux.
solar_flux_lookup(solar_flux, wavelen) Given the observed solar flux at 10.7 cm wavelength, returns the solar flux at the given radar wavelength
scanning_losses(angle_step, beamwidth) Given the antenna beam width and the integration angle, compute the losses due to the fact that the sun is not a point target and the antenna is scanning
smooth_masked(raw_data[, wind_len, …]) smoothes the data using a rolling window.
sun_position_pysolar(dt, lat, lon[, elevation]) obtains the sun position in antenna coordinates using the pysolar library.
sun_position_mfr(dt, lat_deg, lon_deg[, …]) Calculate the sun position for the given time (dt) at the given position (lat, lon).
class pyart.correct.GateFilter(radar, exclude_based=True)[source]

Bases: object

A class for building a boolean arrays for filtering gates based on a set of condition typically based on the values in the radar fields. These filter can be used in various algorithms and calculations within Py-ART.

See pyart.correct.GateFilter.exclude_below() for method parameter details.

Parameters:
radar : Radar

Radar object from which gate filter will be build.

exclude_based : bool, optional

True, the default and suggested method, will begin with all gates included and then use the exclude methods to exclude gates based on conditions. False will begin with all gates excluded from which a set of gates to include should be set using the include methods.

Examples

>>> import pyart
>>> radar = pyart.io.read('radar_file.nc')
>>> gatefilter = pyart.correct.GateFilter(radar)
>>> gatefilter.exclude_below('reflectivity', 10)
>>> gatefilter.exclude_below('normalized_coherent_power', 0.75)
Attributes:
gate_excluded : array, dtype=bool

Return a copy of the excluded gates.

gate_included : array, dtype=bool

Return a copy of the included gates.

Methods

copy(self) Return a copy of the gatefilter.
exclude_above(self, field, value[, …]) Exclude gates where a given field is above a given value.
exclude_all(self) Exclude all gates.
exclude_below(self, field, value[, …]) Exclude gates where a given field is below a given value.
exclude_equal(self, field, value[, …]) Exclude gates where a given field is equal to a value.
exclude_gates(self, mask[, exclude_masked, op]) Exclude gates where a given mask is equal True.
exclude_inside(self, field, v1, v2[, …]) Exclude gates where a given field is inside a given interval.
exclude_invalid(self, field[, …]) Exclude gates where an invalid value occurs in a field (NaNs or infs).
exclude_masked(self, field[, exclude_masked, op]) Exclude gates where a given field is masked.
exclude_none(self) Exclude no gates, include all gates.
exclude_not_equal(self, field, value[, …]) Exclude gates where a given field is not equal to a value.
exclude_outside(self, field, v1, v2[, …]) Exclude gates where a given field is outside a given interval.
exclude_transition(self[, trans_value, …]) Exclude all gates in rays marked as in transition between sweeps.
include_above(self, field, value[, …]) Include gates where a given field is above a given value.
include_all(self) Include all gates.
include_below(self, field, value[, …]) Include gates where a given field is below a given value.
include_equal(self, field, value[, …]) Include gates where a given field is equal to a value.
include_gates(self, mask[, exclude_masked, op]) Include gates where a given mask is equal True.
include_inside(self, field, v1, v2[, …]) Include gates where a given field is inside a given interval.
include_none(self) Include no gates, exclude all gates.
include_not_equal(self, field, value[, …]) Include gates where a given field is not equal to a value.
include_not_masked(self, field[, …]) Include gates where a given field in not masked.
include_not_transition(self[, trans_value, …]) Include all gates in rays not marked as in transition between sweeps.
include_outside(self, field, v1, v2[, …]) Include gates where a given field is outside a given interval.
include_valid(self, field[, exclude_masked, op]) Include gates where a valid value occurs in a field (not NaN or inf).
__class__

alias of builtins.type

__delattr__(self, name, /)

Implement delattr(self, name).

__dict__ = mappingproxy({'__module__': 'pyart.filters.gatefilter', '__doc__': "\n A class for building a boolean arrays for filtering gates based on\n a set of condition typically based on the values in the radar fields.\n These filter can be used in various algorithms and calculations within\n Py-ART.\n\n See :py:func:`pyart.correct.GateFilter.exclude_below` for method\n parameter details.\n\n Parameters\n ----------\n radar : Radar\n Radar object from which gate filter will be build.\n exclude_based : bool, optional\n True, the default and suggested method, will begin with all gates\n included and then use the exclude methods to exclude gates based on\n conditions. False will begin with all gates excluded from which\n a set of gates to include should be set using the include methods.\n\n Attributes\n ----------\n gate_excluded : array, dtype=bool\n Boolean array indicating if a gate should be excluded from a\n calculation. Elements marked True indicate the corresponding gate\n should be excluded. Those marked False should be included.\n This is read-only attribute, any changes to the array will NOT\n be reflected in gate_included and will be lost when the attribute is\n accessed again.\n gate_included : array, dtype=bool\n Boolean array indicating if a gate should be included in a\n calculation. Elements marked True indicate the corresponding gate\n should be include. Those marked False should be excluded.\n This is read-only attribute, any changes to the array will NOT\n be reflected in gate_excluded and will be lost when the attribute is\n accessed again.\n\n Examples\n --------\n >>> import pyart\n >>> radar = pyart.io.read('radar_file.nc')\n >>> gatefilter = pyart.correct.GateFilter(radar)\n >>> gatefilter.exclude_below('reflectivity', 10)\n >>> gatefilter.exclude_below('normalized_coherent_power', 0.75)\n\n ", '__init__': <function GateFilter.__init__>, 'copy': <function GateFilter.copy>, 'gate_included': <property object>, 'gate_excluded': <property object>, '_get_fdata': <function GateFilter._get_fdata>, '_merge': <function GateFilter._merge>, 'exclude_transition': <function GateFilter.exclude_transition>, 'exclude_below': <function GateFilter.exclude_below>, 'exclude_above': <function GateFilter.exclude_above>, 'exclude_inside': <function GateFilter.exclude_inside>, 'exclude_outside': <function GateFilter.exclude_outside>, 'exclude_equal': <function GateFilter.exclude_equal>, 'exclude_not_equal': <function GateFilter.exclude_not_equal>, 'exclude_all': <function GateFilter.exclude_all>, 'exclude_none': <function GateFilter.exclude_none>, 'exclude_masked': <function GateFilter.exclude_masked>, 'exclude_invalid': <function GateFilter.exclude_invalid>, 'exclude_gates': <function GateFilter.exclude_gates>, 'include_not_transition': <function GateFilter.include_not_transition>, 'include_below': <function GateFilter.include_below>, 'include_above': <function GateFilter.include_above>, 'include_inside': <function GateFilter.include_inside>, 'include_outside': <function GateFilter.include_outside>, 'include_equal': <function GateFilter.include_equal>, 'include_not_equal': <function GateFilter.include_not_equal>, 'include_all': <function GateFilter.include_all>, 'include_none': <function GateFilter.include_none>, 'include_not_masked': <function GateFilter.include_not_masked>, 'include_valid': <function GateFilter.include_valid>, 'include_gates': <function GateFilter.include_gates>, '__dict__': <attribute '__dict__' of 'GateFilter' objects>, '__weakref__': <attribute '__weakref__' of 'GateFilter' objects>})
__dir__(self, /)

Default dir() implementation.

__eq__(self, value, /)

Return self==value.

__format__(self, format_spec, /)

Default object formatter.

__ge__(self, value, /)

Return self>=value.

__getattribute__(self, name, /)

Return getattr(self, name).

__gt__(self, value, /)

Return self>value.

__hash__(self, /)

Return hash(self).

__init__(self, radar, exclude_based=True)[source]

initialize

__init_subclass__()

This method is called when a class is subclassed.

The default implementation does nothing. It may be overridden to extend subclasses.

__le__(self, value, /)

Return self<=value.

__lt__(self, value, /)

Return self<value.

__module__ = 'pyart.filters.gatefilter'
__ne__(self, value, /)

Return self!=value.

__new__(*args, **kwargs)

Create and return a new object. See help(type) for accurate signature.

__reduce__(self, /)

Helper for pickle.

__reduce_ex__(self, protocol, /)

Helper for pickle.

__repr__(self, /)

Return repr(self).

__setattr__(self, name, value, /)

Implement setattr(self, name, value).

__sizeof__(self, /)

Size of object in memory, in bytes.

__str__(self, /)

Return str(self).

__subclasshook__()

Abstract classes can override this to customize issubclass().

This is invoked early on by abc.ABCMeta.__subclasscheck__(). It should return True, False or NotImplemented. If it returns NotImplemented, the normal algorithm is used. Otherwise, it overrides the normal algorithm (and the outcome is cached).

__weakref__

list of weak references to the object (if defined)

_get_fdata(self, field)[source]

Check that the field exists and retrieve field data.

_merge(self, marked, op, exclude_masked)[source]

Merge an array of marked gates with the exclude array.

copy(self)[source]

Return a copy of the gatefilter.

exclude_above(self, field, value, exclude_masked=True, op='or', inclusive=False)[source]

Exclude gates where a given field is above a given value.

exclude_all(self)[source]

Exclude all gates.

exclude_below(self, field, value, exclude_masked=True, op='or', inclusive=False)[source]

Exclude gates where a given field is below a given value.

Parameters:
field : str

Name of field compared against the value.

value : float

Gates with a value below this value in the specified field will be marked for exclusion in the filter.

exclude_masked : bool, optional

True to filter masked values in the specified field if the data is a masked array, False to include any masked values.

op : {‘and’, ‘or’, ‘new’}

Operation to perform when merging the existing set of excluded gates with the excluded gates from the current operation. ‘and’ will perform a logical AND operation, ‘or’ a logical OR, and ‘new’ will replace the existing excluded gates with the one generated here. ‘or’, the default for exclude methods, is typically desired when building up a set of conditions for excluding gates where the desired effect is to exclude gates which meet any of the conditions. ‘and’, the default for include methods, is typically desired when building up a set of conditions where the desired effect is to include gates which meet any of the conditions. Note that the ‘and’ method MAY results in including gates which have previously been excluded because they were masked or invalid.

inclusive : bool

Indicates whether the specified value should also be excluded.

exclude_equal(self, field, value, exclude_masked=True, op='or')[source]

Exclude gates where a given field is equal to a value.

exclude_gates(self, mask, exclude_masked=True, op='or')[source]

Exclude gates where a given mask is equal True.

Parameters:
mask : numpy array

Boolean numpy array with same shape as a field array.

exclude_masked : bool, optional

True to filter masked values in the specified mask if it is a masked array, False to include any masked values.

op : {‘and’, ‘or’, ‘new’}

Operation to perform when merging the existing set of excluded gates with the excluded gates from the current operation. ‘and’ will perform a logical AND operation, ‘or’ a logical OR, and ‘new’ will replace the existing excluded gates with the one generated here. ‘or’, the default for exclude methods, is typically desired when building up a set of conditions for excluding gates where the desired effect is to exclude gates which meet any of the conditions. ‘and’, the default for include methods, is typically desired when building up a set of conditions where the desired effect is to include gates which meet any of the conditions. Note that the ‘and’ method MAY results in including gates which have previously been excluded because they were masked or invalid.

exclude_inside(self, field, v1, v2, exclude_masked=True, op='or', inclusive=True)[source]

Exclude gates where a given field is inside a given interval.

exclude_invalid(self, field, exclude_masked=True, op='or')[source]

Exclude gates where an invalid value occurs in a field (NaNs or infs).

exclude_masked(self, field, exclude_masked=True, op='or')[source]

Exclude gates where a given field is masked.

exclude_none(self)[source]

Exclude no gates, include all gates.

exclude_not_equal(self, field, value, exclude_masked=True, op='or')[source]

Exclude gates where a given field is not equal to a value.

exclude_outside(self, field, v1, v2, exclude_masked=True, op='or', inclusive=False)[source]

Exclude gates where a given field is outside a given interval.

exclude_transition(self, trans_value=1, exclude_masked=True, op='or')[source]

Exclude all gates in rays marked as in transition between sweeps.

Exclude all gates in rays marked as “in transition” by the antenna_transition attribute of the radar used to construct the filter. If no antenna transition information is available no gates are excluded.

Parameters:
trans_value : int, optional

Value used in the antenna transition data to indicate that the instrument was between sweeps (in transition) during the collection of a specific ray. Typically a value of 1 is used to indicate this transition and the default can be used in these cases.

exclude_masked : bool, optional

True to filter masked values in antenna_transition if the data is a masked array, False to include any masked values.

op : {‘and’, ‘or’, ‘new’}

Operation to perform when merging the existing set of excluded gates with the excluded gates from the current operation. ‘and’ will perform a logical AND operation, ‘or’ a logical OR, and ‘new’ will replace the existing excluded gates with the one generated here. ‘or’, the default for exclude methods, is typically desired when building up a set of conditions for excluding gates where the desired effect is to exclude gates which meet any of the conditions. ‘and’, the default for include methods, is typically desired when building up a set of conditions where the desired effect is to include gates which meet any of the conditions. Note that the ‘and’ method MAY results in including gates which have previously been excluded because they were masked or invalid.

gate_excluded

Return a copy of the excluded gates.

gate_included

Return a copy of the included gates.

include_above(self, field, value, exclude_masked=True, op='and', inclusive=False)[source]

Include gates where a given field is above a given value.

include_all(self)[source]

Include all gates.

include_below(self, field, value, exclude_masked=True, op='and', inclusive=False)[source]

Include gates where a given field is below a given value.

include_equal(self, field, value, exclude_masked=True, op='and')[source]

Include gates where a given field is equal to a value.

include_gates(self, mask, exclude_masked=True, op='and')[source]

Include gates where a given mask is equal True.

Parameters:
mask : numpy array

Boolean numpy array with same shape as a field array.

exclude_masked : bool, optional

True to filter masked values in the specified mask if it is a masked array, False to include any masked values.

op : {‘and’, ‘or’, ‘new’}

Operation to perform when merging the existing set of excluded gates with the excluded gates from the current operation. ‘and’ will perform a logical AND operation, ‘or’ a logical OR, and ‘new’ will replace the existing excluded gates with the one generated here. ‘or’, the default for exclude methods, is typically desired when building up a set of conditions for excluding gates where the desired effect is to exclude gates which meet any of the conditions. ‘and’, the default for include methods, is typically desired when building up a set of conditions where the desired effect is to include gates which meet any of the conditions. Note that the ‘or’ method MAY results in excluding gates which have previously been included.

include_inside(self, field, v1, v2, exclude_masked=True, op='and', inclusive=True)[source]

Include gates where a given field is inside a given interval.

include_none(self)[source]

Include no gates, exclude all gates.

include_not_equal(self, field, value, exclude_masked=True, op='and')[source]

Include gates where a given field is not equal to a value.

include_not_masked(self, field, exclude_masked=True, op='and')[source]

Include gates where a given field in not masked.

include_not_transition(self, trans_value=0, exclude_masked=True, op='and')[source]

Include all gates in rays not marked as in transition between sweeps.

Include all gates in rays not marked as “in transition” by the antenna_transition attribute of the radar used to construct the filter. If no antenna transition information is available all gates are included.

Parameters:
trans_value : int, optional

Value used in the antenna transition data to indicate that the instrument is not between sweeps (in transition) during the collection of a specific ray. Typically a value of 0 is used to indicate no transition and the default can be used in these cases.

exclude_masked : bool, optional

True to filter masked values in antenna_transition if the data is a masked array, False to include any masked values.

op : {‘and’, ‘or’, ‘new’}

Operation to perform when merging the existing set of excluded gates with the excluded gates from the current operation. ‘and’ will perform a logical AND operation, ‘or’ a logical OR, and ‘new’ will replace the existing excluded gates with the one generated here. ‘or’, the default for exclude methods, is typically desired when building up a set of conditions for excluding gates where the desired effect is to exclude gates which meet any of the conditions. ‘and’, the default for include methods, is typically desired when building up a set of conditions where the desired effect is to include gates which meet any of the conditions. Note that the ‘or’ method MAY results in excluding gates which have previously been included.

include_outside(self, field, v1, v2, exclude_masked=True, op='and', inclusive=False)[source]

Include gates where a given field is outside a given interval.

include_valid(self, field, exclude_masked=True, op='and')[source]

Include gates where a valid value occurs in a field (not NaN or inf).

pyart.correct.calculate_attenuation(radar, z_offset, debug=False, doc=15, fzl=4000.0, rhv_min=0.8, ncp_min=0.5, a_coef=0.06, beta=0.8, refl_field=None, ncp_field=None, rhv_field=None, phidp_field=None, spec_at_field=None, corr_refl_field=None)[source]

Calculate the attenuation from a polarimetric radar using Z-PHI method. Parameters ———- radar : Radar

Radar object to use for attenuation calculations. Must have copol_coeff, norm_coherent_power, proc_dp_phase_shift, reflectivity_horizontal fields.
z_offset : float
Horizontal reflectivity offset in dBZ.
debug : bool, optional
True to print debugging information, False supressed this printing.
doc : float, optional
Number of gates at the end of each ray to to remove from the calculation.
fzl : float, optional
Freezing layer, gates above this point are not included in the correction.
rhv_min : float, optional
Minimum copol_coeff value to consider valid.
ncp_min : float, optional
Minimum norm_coherent_power to consider valid.
a_coef : float, optional
A coefficient in attenuation calculation.
beta : float, optional
Beta parameter in attenuation calculation.
refl_field : str, optional
Name of the reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.
phidp_field : str, optional
Name of the differential phase field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.
ncp_field : str, optional
Name of the normalized coherent power field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.
zdr_field : str, optional
Name of the differential reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file. This will only be used if it is available.
spec_at_field : str, optional
Name of the specific attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.
corr_refl_field : str, optional
Name of the corrected reflectivity field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.
Returns:
spec_at : dict

Field dictionary containing the specific attenuation.

cor_z : dict

Field dictionary containing the corrected reflectivity.

References

Gu et al. Polarimetric Attenuation Correction in Heavy Rain at C Band, JAMC, 2011, 50, 39-58.

pyart.correct.calculate_attenuation_philinear(radar, doc=None, fzl=None, pia_coef=None, gatefilter=None, pida_coef=None, refl_field=None, phidp_field=None, zdr_field=None, temp_field=None, iso0_field=None, spec_at_field=None, pia_field=None, corr_refl_field=None, spec_diff_at_field=None, pida_field=None, corr_zdr_field=None, temp_ref='temperature')[source]

Calculate the attenuation and the differential attenuation from a polarimetric radar using linear dependece with PhiDP. The attenuation is computed up to a user defined freezing level height, where temperatures in a temperature field are positive or where the height relative to the iso0 is 0. The coefficients are either user-defined or radar frequency dependent.

Parameters:
radar : Radar

Radar object to use for attenuation calculations. Must have phidp and refl fields.

doc : float, optional

Number of gates at the end of each ray to to remove from the calculation.

fzl : float, optional

Freezing layer, gates above this point are not included in the correction.

gatefilter : GateFilter, optional

The gates to exclude from the calculation. This, combined with the gates above fzl, will be excluded from the correction. Set to None to not use a gatefilter.

pia_coef : float, optional

Coefficient in path integrated attenuation calculation

pida_coeff : float, optional

Coefficient in path integrated differential attenuation calculation

refl_field : str, optional

Name of the reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

phidp_field : str, optional

Name of the differential phase field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

zdr_field : str, optional

Name of the differential reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file. This will only be used if it is available.

temp_field : str, optional

Name of the temperature field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

iso0_field : str, optional

Name of the field for the height above the 0C isotherm for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file. This will only be used if it is available.

spec_at_field : str, optional

Name of the specific attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

pia_field : str, optional

Name of the path integrated attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

corr_refl_field : str, optional

Name of the corrected reflectivity field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

spec_diff_at_field : str, optional

Name of the specific differential attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file. This will only be calculated if ZDR is available.

corr_zdr_field : str, optional

Name of the corrected differential reflectivity field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file. This will only be calculated if ZDR is available.

temp_ref : str, optional

The field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl.

Returns:
spec_at : dict

Field dictionary containing the specific attenuation.

pia_dict : dict

Field dictionary containing the path integrated attenuation.

cor_z : dict

Field dictionary containing the corrected reflectivity.

spec_diff_at : dict

Field dictionary containing the specific differential attenuation.

pida_dict : dict

Field dictionary containing the path integrated differential attenuation.

cor_zdr : dict

Field dictionary containing the corrected differential reflectivity.

pyart.correct.calculate_attenuation_zphi(radar, doc=None, fzl=None, smooth_window_len=5, gatefilter=None, a_coef=None, beta=None, c=None, d=None, refl_field=None, phidp_field=None, zdr_field=None, temp_field=None, iso0_field=None, spec_at_field=None, pia_field=None, corr_refl_field=None, spec_diff_at_field=None, pida_field=None, corr_zdr_field=None, temp_ref='temperature')[source]

Calculate the attenuation and the differential attenuation from a polarimetric radar using Z-PHI method.. The attenuation is computed up to a user defined freezing level height or up to where temperatures in a temperature field are positive. The coefficients are either user-defined or radar frequency dependent.

Parameters:
radar : Radar

Radar object to use for attenuation calculations. Must have phidp and refl fields.

doc : float, optional

Number of gates at the end of each ray to to remove from the calculation.

fzl : float, optional

Freezing layer, gates above this point are not included in the correction.

gatefilter : GateFilter, optional

The gates to exclude from the calculation. This, combined with the gates above fzl, will be excluded from the correction. Set to None to not use a gatefilter.

smooth_window_len : int, optional

Size, in range bins, of the smoothing window

a_coef : float, optional

A coefficient in attenuation calculation.

beta : float, optional

Beta parameter in attenuation calculation.

c, d : float, optional

coefficient and exponent of the power law that relates attenuation with differential attenuation

refl_field : str, optional

Name of the reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

phidp_field : str, optional

Name of the differential phase field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

zdr_field : str, optional

Name of the differential reflectivity field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file. This will only be used if it is available.

temp_field : str, optional

Name of the temperature field used for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

iso0_field : str, optional

Name of the field for the height above the 0C isotherm for the attenuation correction. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file. This will only be used if it is available.

spec_at_field : str, optional

Name of the specific attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

pia_field : str, optional

Name of the path integrated attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

corr_refl_field : str, optional

Name of the corrected reflectivity field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file.

spec_diff_at_field : str, optional

Name of the specific differential attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file. This will only be calculated if ZDR is available.

pida_field : str, optional

Name of the path integrated differential attenuation field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file. This will only be calculated if ZDR is available.

corr_zdr_field : str, optional

Name of the corrected differential reflectivity field that will be used to fill in the metadata for the returned fields. A value of None for any of these parameters will use the default field names as defined in the Py-ART configuration file. This will only be calculated if ZDR is available.

temp_ref : str, optional

the field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl

Returns:
spec_at : dict

Field dictionary containing the specific attenuation.

pia_dict : dict

Field dictionary containing the path integrated attenuation.

cor_z : dict

Field dictionary containing the corrected reflectivity.

spec_diff_at : dict

Field dictionary containing the specific differential attenuation.

pida_dict : dict

Field dictionary containing the path integrated differential attenuation.

cor_zdr : dict

Field dictionary containing the corrected differential reflectivity.

References

Gu et al. Polarimetric Attenuation Correction in Heavy Rain at C Band, JAMC, 2011, 50, 39-58.

Ryzhkov et al. Potential Utilization of Specific Attenuation for Rainfall Estimation, Mitigation of Partial Beam Blockage, and Radar Networking, JAOT, 2014, 31, 599-619.

pyart.correct.correct_bias(radar, bias=0.0, field_name=None)[source]

Corrects a radar data bias. If field name is none the correction is applied to horizontal reflectivity by default.

Parameters:
radar : Radar

Radar object.

bias : float, optional

The bias magnitude.

field_name: str, optional

Names of the field to be corrected.

Returns:
corrected_field : dict

The corrected field

pyart.correct.correct_noise_rhohv(radar, urhohv_field=None, snr_field=None, zdr_field=None, nh_field=None, nv_field=None, rhohv_field=None)[source]

Corrects RhoHV for noise according to eq. 6 in Gourley et al. 2006. This correction should only be performed if noise has not been subtracted from the signal during the moments computation.

Parameters:
radar : Radar

Radar object.

urhohv_field : str, optional

Name of the RhoHV uncorrected for noise field.

snr_field, zdr_field, nh_field, nv_field : str, optional

Names of the SNR, ZDR, horizontal channel noise in dBZ and vertical channel noise in dBZ used to correct RhoHV.

rhohv_field : str, optional

Name of the rhohv field to output.

Returns:
rhohv : dict

Noise corrected RhoHV field.

References

Gourley et al. Data Quality of the Meteo-France C-Band Polarimetric Radar, JAOT, 23, 1340-1356

pyart.correct.correct_sys_phase(radar, ind_rmin=10, ind_rmax=500, min_rcons=11, zmin=20.0, zmax=40.0, psidp_field=None, refl_field=None, phidp_field=None)[source]

correction of the system offset. Public method

Parameters:
radar : Radar

Radar object for which to determine the system phase.

ind_rmin, ind_rmax : int

Min and max range index where to look for continuous precipitation

min_rcons : int

The minimum number of consecutive gates to consider it a rain cell.

zmin, zmax : float

Minimum and maximum reflectivity to consider it a rain cell

psidp_field : str

Field name within the radar object which represent the differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

refl_field : str

Field name within the radar object which represent the reflectivity. A value of None will use the default field name as defined in the Py-ART configuration file.

phidp_field : str

Field name within the radar object which represent the corrected differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

Returns:
phidp_dict : dict

The corrected phidp field

pyart.correct.correct_visibility(radar, vis_field=None, field_name=None)[source]

Corrects the reflectivity according to visibility. Applied to horizontal reflectivity by default

Parameters:
radar : Radar

radar object

vis_field : str

the name of the visibility field

field_name: str

names of the field to be corrected

Returns:
corrected_field : dict

The corrected field

pyart.correct.dealias_fourdd(radar, last_radar=None, sonde_profile=None, gatefilter=False, filt=1, rsl_badval=131072.0, keep_original=False, set_limits=True, vel_field=None, corr_vel_field=None, last_vel_field=None, debug=False, max_shear=0.05, sign=1, **kwargs)[source]

Dealias Doppler velocities using the 4DD algorithm.

Dealias the Doppler velocities field using the University of Washington 4DD algorithm utilizing information from a previous volume scan and/or sounding data. Either last_radar or sonde_profile must be provided. For best results provide both a previous volume scan and sounding data. Radar and last_radar must contain the same number of rays per sweep.

Additional arguments are passed to _fourdd_interface.fourdd_dealias(). These can be used to fine tune the behavior of the FourDD algorithm. See the documentation of Other Parameters for details. For the default values of these parameters see the documentation of _fourdd_interface.fourdd_dealias().

Parameters:
radar : Radar

Radar object to use for dealiasing. Must have a Nyquist defined in the instrument_parameters attribute and have a reflectivity_horizontal and mean_doppler_velocity fields.

last_radar : Radar, optional

The previous radar volume, which has been successfully dealiased. Using a previous volume as an initial condition can greatly improve the dealiasing, and represents the final dimension in the 4DD algorithm.

sonde_profile : HorizontalWindProfile, optional

Profile of horizontal winds from a sonding used for the initial condition of the dealiasing.

Returns:
vr_corr : dict

Field dictionary containing dealiased Doppler velocities. Dealiased array is stored under the ‘data’ key.

Other Parameters:
 
gatefilter : GateFilter, optional.

A GateFilter instance which specifies which gates should be ignored when performing velocity dealiasing. A value of None will create this filter from the radar moments using any additional arguments by passing them to moment_based_gate_filter(). The default value assumes all gates are valid.

filt : int, optional

Flag controlling Bergen and Albers filter, 1 = yes, 0 = no.

rsl_badval : float, optional

Value which represents a bad value in RSL.

keep_original : bool, optional

True to keep original doppler velocity values when the dealiasing procedure fails, otherwise these gates will be masked. NaN values are still masked.

set_limits : bool, optional

True to set valid_min and valid_max elements in the returned dictionary. False will not set these dictionary elements.

vel_field : str, optional

Field in radar to use as the Doppler velocities during dealiasing. None will use the default field name from the Py-ART configuration file.

corr_vel_field : str, optional

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

last_vel_field : str, optional

Name to use for the dealiased Doppler velocity field metadata in last_radar. None will use the corr_vel_field name.

maxshear : float, optional

Maximum vertical shear which will be incorporated into the created volume from the sounding data. Parameter not used when no sounding data is provided.

sign : int, optional

Sign convention which the radial velocities in the volume created from the sounding data will will. This should match the convention used in the radar data. A value of 1 represents when positive values velocities are towards the radar, -1 represents when negative velocities are towards the radar.

compthresh : float, optional

Fraction of the Nyquist velocity to use as a threshold when performing continuity (initial) dealiasing. Velocities differences above this threshold will not be marked as gate from which to begin unfolding during spatial dealiasing.

compthresh2 : float, optional

The same as compthresh but the value used during the second pass of dealiasing. This second pass is only performed in both a sounding and last volume are provided.

thresh : float, optional

Fraction of the Nyquist velocity to use as a threshold when performing spatial dealiasing. Horizontally adjacent gates with velocities above this threshold will count against assigning the gate in question the velocity value being tested.

ckval : float, optional

When the absolute value of the velocities are below this value they will not be marked as gates from which to begin unfolding during spatial dealiasing.

stdthresh : float, optional

Fraction of the Nyquist velocity to use as a standard deviation threshold in the window dealiasing portion of the algorithm.

epsilon : float, optional

Difference used when comparing a value to missing value, changing this from the default is not recommended.

maxcount : int, optional

Maximum allowed number of fold allowed when unfolding velocities.

pass2 : int, optional

Controls weather unfolded gates should be removed (a value of 0) or retained for unfolding during the second pass (a value of 1) when both a sounding volume and last volume are provided.

rm : int, optional

Determines what should be done with gates that are left unfolded after the first pass of dealiasing. A value of 1 will remove these gates, a value of 0 sets these gates to their initial velocity. If both a sounding volume and last volume are provided this parameter is ignored.

proximity : int, optional

Number of gates and rays to include of either side of the current gate during window dealiasing. This value may be doubled in cases where a standard sized window does not capture a sufficient number of good valued gates.

mingood : int, optional

Number of good valued gates required within the window before the current gate will be unfolded.

ba_mincount : int, optional

Number of neighbors required during Bergen and Albers filter for a given gate to be included, must be between 1 and 8, 5 recommended.

ba_edgecount : int, optional

Same as ba_mincount but used at ray edges, must be between 1 and 5, 3 recommended.

debug : bool, optional

Set True to return RSL Volume objects for debugging: usuccess, radialVelVolume, lastVelVolume, unfoldedVolume, sondVolume

Notes

Due to limitations in the C code do not call with sounding arrays over 999 elements long.

References

C. N. James and R. A Houze Jr, A Real-Time Four-Dimensional Doppler Dealising Scheme, Journal of Atmospheric and Oceanic Technology, 2001, 18, 1674.

pyart.correct.dealias_region_based(radar, ref_vel_field=None, interval_splits=3, interval_limits=None, skip_between_rays=100, skip_along_ray=100, centered=True, nyquist_vel=None, check_nyquist_uniform=True, gatefilter=False, rays_wrap_around=None, keep_original=False, set_limits=True, vel_field=None, corr_vel_field=None, **kwargs)[source]

Dealias Doppler velocities using a region based algorithm.

Performs Doppler velocity dealiasing by finding regions of similar velocities and unfolding and merging pairs of regions until all regions are unfolded. Unfolding and merging regions is accomplished by modeling the problem as a dynamic network reduction.

Parameters:
radar : Radar

Radar object containing Doppler velocities to dealias.

ref_vel_field : str or None, optional

Field in radar containing a reference velocity field used to anchor the unfolded velocities once the algorithm completes. Typically this field is created by simulating the radial velocities from wind data from an atmospheric sonding using pyart.util.simulated_vel_from_profile().

interval_splits : int, optional

Number of segments to split the nyquist interval into when finding regions of similar velocity. More splits creates a larger number of initial regions which takes longer to process but may result in better dealiasing. The default value of 3 seems to be a good compromise between performance and artifact free dealiasing. This value is not used if the interval_limits parameter is not None.

interval_limits : array like or None, optional

Velocity limits used for finding regions of similar velocity. Should cover the entire nyquist interval. None, the default value, will split the Nyquist interval into interval_splits equal sized intervals.

skip_between_rays, skip_along_ray : int, optional

Maximum number of filtered gates to skip over when joining regions, gaps between region larger than this will not be connected. Parameters specify the maximum number of filtered gates between and along a ray. Set these parameters to 0 to disable unfolding across filtered gates.

centered : bool, optional

True to apply centering to each sweep after the dealiasing algorithm so that the average number of unfolding is near 0. False does not apply centering which may results in individual sweeps under or over folded by the nyquist interval.

nyquist_velocity : array like or float, optional

Nyquist velocity in unit identical to those stored in the radar’s velocity field, either for each sweep or a single value which will be used for all sweeps. None will attempt to determine this value from the Radar object.

check_nyquist_uniform : bool, optional

True to check if the Nyquist velocities are uniform for all rays within a sweep, False will skip this check. This parameter is ignored when the nyquist_velocity parameter is not None.

gatefilter : GateFilter, None or False, optional.

A GateFilter instance which specified which gates should be ignored when performing de-aliasing. A value of None created this filter from the radar moments using any additional arguments by passing them to moment_based_gate_filter(). False, the default, disables filtering including all gates in the dealiasing.

rays_wrap_around : bool or None, optional

True when the rays at the beginning of the sweep and end of the sweep should be interpreted as connected when de-aliasing (PPI scans). False if they edges should not be interpreted as connected (other scan types). None will determine the correct value from the radar scan type.

keep_original : bool, optional

True to retain the original Doppler velocity values at gates where the dealiasing procedure fails or was not applied. False does not replacement and these gates will be masked in the corrected velocity field.

set_limits : bool, optional

True to set valid_min and valid_max elements in the returned dictionary. False will not set these dictionary elements.

vel_field : str, optional

Field in radar to use as the Doppler velocities during dealiasing. None will use the default field name from the Py-ART configuration file.

corr_vel_field : str, optional

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

Returns:
corr_vel : dict

Field dictionary containing dealiased Doppler velocities. Dealiased array is stored under the ‘data’ key.

pyart.correct.dealias_unwrap_phase(radar, unwrap_unit='sweep', nyquist_vel=None, check_nyquist_uniform=True, gatefilter=False, rays_wrap_around=None, keep_original=False, set_limits=True, vel_field=None, corr_vel_field=None, skip_checks=False, **kwargs)[source]

Dealias Doppler velocities using multi-dimensional phase unwrapping.

Parameters:
radar : Radar

Radar object containing Doppler velocities to dealias.

unwrap_unit : {‘ray’, ‘sweep’, ‘volume’}, optional

Unit to unwrap independently. ‘ray’ will unwrap each ray individually, ‘sweep’ each sweep, and ‘volume’ will unwrap the entire volume in a single pass. ‘sweep’, the default, often gives superior results when the lower sweeps of the radar volume are contaminated by clutter. ‘ray’ does not use the gatefilter parameter and rays where gates ared masked will result in poor dealiasing for that ray.

nyquist_velocity : array like or float, optional

Nyquist velocity in unit identical to those stored in the radar’s velocity field, either for each sweep or a single value which will be used for all sweeps. None will attempt to determine this value from the Radar object. The Nyquist velocity of the first sweep is used for all dealiasing unless the unwrap_unit is ‘sweep’ when the velocities of each sweep are used.

check_nyquist_uniform : bool, optional

True to check if the Nyquist velocities are uniform for all rays within a sweep, False will skip this check. This parameter is ignored when the nyquist_velocity parameter is not None.

gatefilter : GateFilter, None or False, optional.

A GateFilter instance which specified which gates should be ignored when performing de-aliasing. A value of None created this filter from the radar moments using any additional arguments by passing them to moment_based_gate_filter(). False, the default, disables filtering including all gates in the dealiasing.

rays_wrap_around : bool or None, optional

True when the rays at the beginning of the sweep and end of the sweep should be interpreted as connected when de-aliasing (PPI scans). False if they edges should not be interpreted as connected (other scan types). None will determine the correct value from the radar scan type.

keep_original : bool, optional

True to retain the original Doppler velocity values at gates where the dealiasing procedure fails or was not applied. False does not replacement and these gates will be masked in the corrected velocity field.

set_limits : bool, optional

True to set valid_min and valid_max elements in the returned dictionary. False will not set these dictionary elements.

vel_field : str, optional

Field in radar to use as the Doppler velocities during dealiasing. None will use the default field name from the Py-ART configuration file.

corr_vel_field : str, optional

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

skip_checks : bool

True to skip checks verifing that an appropiate unwrap_unit is selected, False retains these checked. Setting this parameter to True is not recommended and is only offered as an option for extreme cases.

Returns:
corr_vel : dict

Field dictionary containing dealiased Doppler velocities. Dealiased array is stored under the ‘data’ key.

References

[1]Miguel Arevallilo Herraez, David R. Burton, Michael J. Lalor, and Munther A. Gdeisat, “Fast two-dimensional phase-unwrapping algorithm based on sorting by reliability following a noncontinuous path”, Journal Applied Optics, Vol. 41, No. 35 (2002) 7437,
[2]Abdul-Rahman, H., Gdeisat, M., Burton, D., & Lalor, M., “Fast three-dimensional phase-unwrapping algorithm based on sorting by reliability following a non-continuous path. In W. Osten, C. Gorecki, & E. L. Novak (Eds.), Optical Metrology (2005) 32–40, International Society for Optics and Photonics.
pyart.correct.despeckle_field(radar, field, label_dict=None, threshold=-100, size=10, gatefilter=None, delta=5.0)[source]

Despeckle a radar volume by identifying small objects in each scan and masking them out. User can define which field to investigate, as well as various thresholds to use on that field and any objects found within. Requires scipy to be installed, and returns a GateFilter object.

Parameters:
radar : pyart.core.Radar object

Radar object to query.

field : str

Name of field to investigate for speckles.

Returns:
gatefilter : pyart.filters.GateFilter object

Py-ART GateFilter object that includes the despeckling mask

Other Parameters:
 
label_dict : dict or None, optional

Dictionary that is produced by find_objects. If None, find_objects will be called to produce it.

threshold : int or float, or 2-element tuple of ints or floats

Threshold values above (if single value) or between (if tuple) for objects to be identified. Default value assumes reflectivity.

size : int, optional

Number of contiguous gates in an object, below which it is a speckle.

gatefilter : None or pyart.filters.GateFilter object

Py-ART GateFilter object to which to add the despeckling mask. The GateFilter object will be permanently modified with the new filtering. If None, creates a new GateFilter.

delta : int or float, optional

Size of allowable gap near PPI edges, in deg, to consider it full 360. If gap is small, then PPI edges will be checked for matching objects.

pyart.correct.det_sys_phase_ray(radar, ind_rmin=10, ind_rmax=500, min_rcons=11, zmin=20.0, zmax=40.0, phidp_field=None, refl_field=None)[source]

Public method Alternative determination of the system phase. Assumes that the valid gates of phidp are only precipitation. A system phase value is found for each ray.

Parameters:
radar : Radar

Radar object for which to determine the system phase.

ind_rmin, ind_rmax : int

Min and max range index where to look for continuous precipitation

min_rcons : int

The minimum number of consecutive gates to consider it a rain cell.

zmin, zmax : float

The minimum and maximum reflectivity to consider the radar bin suitable precipitation

phidp_field : str

Field name within the radar object which represent the differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

refl_field : str

Field name within the radar object which represent the reflectivity. A value of None will use the default field name as defined in the Py-ART configuration file.

Returns:
phidp0_dict : dict

Estimate of the system phase at each ray and metadata

first_gates_dict : dict

The first gate where PhiDP is valid and metadata

pyart.correct.est_rhohv_rain(radar, ind_rmin=10, ind_rmax=500, zmin=20.0, zmax=40.0, thickness=700.0, doc=None, fzl=None, rhohv_field=None, temp_field=None, iso0_field=None, refl_field=None, temp_ref='temperature')[source]

Estimates the quantiles of RhoHV in rain for each sweep

Parameters:
radar : Radar

radar object

ind_rmin, ind_rmax : int

Min and max range index where to look for rain

zmin, zmax : float

The minimum and maximum reflectivity to consider the radar bin suitable rain

thickness : float

Assumed thickness of the melting layer

doc : float

Number of gates at the end of each ray to to remove from the calculation.

fzl : float

Freezing layer, gates above this point are not included in the correction.

temp_field, iso0_field, rhohv_field, refl_field : str

Field names within the radar object which represent the temperature, the height over the iso0, co-polar correlation and reflectivity fields. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_ref : str

the field use as reference for temperature. Can be either temperature or height_over_iso0

Returns:
rhohv_rain_dict : dict

The estimated RhoHV in rain for each sweep and metadata

pyart.correct.est_zdr_precip(radar, ind_rmin=10, ind_rmax=500, zmin=20.0, zmax=22.0, rhohvmin=0.97, phidpmax=10.0, elmax=None, thickness=700.0, doc=None, fzl=None, zdr_field=None, rhohv_field=None, phidp_field=None, temp_field=None, iso0_field=None, refl_field=None, temp_ref='temperature')[source]

Filters out all undesired data to be able to estimate ZDR bias, either in moderate rain or from vertically pointing scans

Parameters:
radar : Radar

radar object

ind_rmin, ind_rmax : int

Min and max range index where to look for rain

zmin, zmax : float

The minimum and maximum reflectivity to consider the radar bin suitable rain

rhohvmin : float

Minimum RhoHV to consider the radar bin suitable rain

phidpmax : float

Maximum PhiDP to consider the radar bin suitable rain

elmax : float

Maximum elevation

thickness : float

Assumed thickness of the melting layer

doc : float

Number of gates at the end of each ray to to remove from the calculation.

fzl : float

Freezing layer, gates above this point are not included in the correction.

zdr_field, rhohv_field, refl_field, phidp_field, temp_field,

iso0_field : str Field names within the radar object which represent the differential reflectivity, co-polar correlation, reflectivity, differential phase, temperature and height relative to the iso0 fields. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_ref : str

the field use as reference for temperature. Can be either temperature, height_over_iso0, fixed_fzl or None

Returns:
zdr_prec_dict : dict

The ZDR data complying with specifications and metadata

pyart.correct.est_zdr_snow(radar, ind_rmin=10, ind_rmax=500, zmin=0.0, zmax=30.0, snrmin=10.0, snrmax=50.0, rhohvmin=0.97, kept_values=[2], phidpmax=10.0, kdpmax=None, tempmin=None, tempmax=None, elmax=None, zdr_field=None, rhohv_field=None, phidp_field=None, temp_field=None, snr_field=None, hydro_field=None, kdp_field=None, refl_field=None)[source]

Filters out all undesired data to be able to estimate ZDR bias in snow

Parameters:
radar : Radar

radar object

ind_rmin, ind_rmax : int

Min and max range index where to look for snow

zmin, zmax : float

The minimum and maximum reflectivity to consider the radar bin suitable snow

snrmin, snrmax : float

The minimum and maximum SNR to consider the radar bin suitable snow

rhohvmin : float

Minimum RhoHV to consider the radar bin suitable snow

kept_values : list of int

The hydrometeor classification values to keep

phidpmax : float

Maximum PhiDP to consider the radar bin suitable snow

kdpmax : float or None

Maximum KDP. If not none this is the maximum KDP value to consider the radar bin suitable snow

tempmin, tempmax : float or None

If not None, the minimum and maximum temperature to consider the radar bin suitable snow

elmax : float

Maximum elevation

zdr_field, rhohv_field, refl_field, phidp_field, kdp_field, temp_field,
snr_field, hydro_field : str

Field names within the radar object which represent the differential reflectivity, co-polar correlation, reflectivity, differential phase, specific differnetial phase, signal to noise ratio, hydrometeor classification and temperature fields. A value of None will use the default field name as defined in the Py-ART configuration file.

Returns:
zdr_snow_dict : dict

The ZDR data complying with specifications and metadata

pyart.correct.find_objects(radar, field, threshold, sweeps=None, smooth=None, gatefilter=None, delta=5.0)[source]

Find objects (i.e., contiguous gates) in one or more sweeps that match thresholds. Filtering & smoothing are available prior to labeling objects. In addition, periodic boundaries are accounted for if they exist (e.g., 360-deg PPIs). Requires scipy to be installed.

Parameters:
radar : pyart.core.Radar object

Radar object to query.

field : str

Name of field to investigate for objects.

threshold : int or float, or 2-element tuple of ints or floats

Threshold values above (if single value) or between (if tuple) for objects to be identified.

Returns:
label_dict : dict

Dictionary that contains all the labeled objects. If this function is performed on the full Radar object, then the dict is ready to be added as a field.

Other Parameters:
 
sweeps : int or array of ints or None, optional

Sweep numbers to examine. If None, all sweeps are examined.

smooth : int or None, optional

Number of gates included in a smoothing box filter along a ray. If None, no smoothing is done prior to labeling objects.

gatefilter : None or pyart.filters.GateFilter object, optional

Py-ART GateFilter object to apply before labeling objects. If None, no filtering will be performed. Note: Filtering always occurs before smoothing.

delta : int or float, optional

Size of allowable gap near PPI edges, in deg, to consider it full 360. If gap is small, then PPI edges will be checked for matching objects along the periodic boundary.

pyart.correct.get_mask_fzl(radar, fzl=None, doc=None, min_temp=0.0, max_h_iso0=0.0, thickness=None, beamwidth=None, temp_field=None, iso0_field=None, temp_ref='temperature')[source]

Constructs a mask to mask data placed thickness m below data at min_temp and beyond.

Parameters:
radar : Radar

The radar object.

fzl : float, optional

Freezing layer, gates above this point are not included in the correction.

doc : float, optional

Number of gates at the end of each ray to to remove from the calculation.

min_temp : float, optional

Minimum temperature below which the data is mask in degrees.

max_h_iso0 : float, optional

Maximum height relative to the iso0 below which the data is mask in meters.

thickness : float, optional

Extent of the layer below the first gate where min_temp is reached that is going to be masked.

beamwidth : float, optional

The radar antenna 3 dB beamwidth.

temp_field: str, optional

The temperature field. A value of None will use the default field name as defined in the Py-ART configuration file. It is going to be used only if available.

iso0_field: str, optional

The field containing the height over the 0C isotherm. A value of None will use the default field name as defined in the Py-ART configuration file. It is going to be used only if available.

temp_ref : str, optional

The field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl.

Returns:
mask_fzl : 2D array

The values that should be masked.

end_gate_arr : 1D array

The index of the last valid gate in the ray.

pyart.correct.get_sun_hits(radar, delev_max=2.0, dazim_max=2.0, elmin=1.0, rmin=50000.0, hmin=10000.0, nbins_min=20, attg=None, max_std_pwr=1.0, max_std_zdr=1.5, sun_position='MF', pwrh_field=None, pwrv_field=None, zdr_field=None)[source]

get data from suspected sun hits. The sun hits are detected using the Hildebrand and Sekhon noise estimation

Parameters:
radar : Radar

radar object

delev_max, dazim_max : float

maximum difference in elevation and azimuth between sun position and antenna pointing

elmin : float

minimum radar elevation angle

rmin : float

minimum range from which we can look for noise [m]

hmin : float

minimum altitude from which we can look for noise [m]. The actual range min will be the minimum between rmin and the range bin higher than hmin.

nbins_min : int

Minimum number of bins with valid data to consider a ray as potentially sun hit

attg : float

gas attenuation coefficient (1-way)

max_std_pwr : float

Maximum standard deviation of the estimated sun power to consider the sun signal valid [dB]

max_std_zdr : float

Maximum standard deviation of the estimated sun ZDR to consider the sun signal valid [dB]

sun_position : str

The function to use to compute the sun position. Can be ‘MF’ or ‘pysolar’

pwrh_field, pwrv_field, zdr_field : str

names of the signal power in dBm for the H and V polarizations and the differential reflectivity

Returns:
sun_hits : dict

a dictionary containing information of the sun hits

new_radar : radar object

radar object containing sweeps that contain sun hits

pyart.correct.get_sun_hits_ivic(radar, delev_max=2.0, dazim_max=2.0, elmin=1.0, npulses_ray=30, flat_reg_wlen=96, nbins_min=800, iterations=10, attg=None, sun_position='MF', max_std_zdr=1.5, pwrh_field=None, pwrv_field=None, zdr_field=None)[source]

get data from suspected sun hits. The sun hits are detected using the Ivic et al. (2003) noise estimator

Parameters:
radar : Radar

radar object

delev_max, dazim_max : float

maximum difference in elevation and azimuth between sun position and antenna pointing

elmin : float

minimum radar elevation angle

npulses_ray : int

Default number of pulses used in the computation of the ray. If the number of pulses is not in radar.instrument_parameters this will be used instead

flat_reg_wlen : int

number of gates considered to find flat regions. The number represents 8 km length with a 83.3 m resolution

nbins_min: int

minimum number of gates with noise to consider the retrieval valid

iterations: int

number of iterations in step 7 of Ivic algorithm

attg : float

gas attenuation coefficient (1-way)

sun_position : str

The function to use to compute the sun position. Can be ‘MF’ or ‘pysolar’

max_std_zdr : float

Maximum standard deviation of the estimated sun ZDR to consider the sun signal valid [dB]

pwrh_field, pwrv_field, zdr_field : str

names of the signal power in dBm for the H and V polarizations and the differential reflectivity

Returns:
sun_hits : dict

a dictionary containing information of the sun hits

new_radar : radar object

radar object containing sweeps that contain sun hits

pyart.correct.moment_based_gate_filter(radar, ncp_field=None, rhv_field=None, refl_field=None, min_ncp=0.5, min_rhv=None, min_refl=-20.0, max_refl=100.0)[source]

Create a filter which removes undesired gates based on moments.

Creates a gate filter in which the following gates are excluded:

  • Gates where the instrument is transitioning between sweeps.
  • Gates where the reflectivity is outside the interval min_refl, max_refl.
  • Gates where the normalized coherent power is below min_ncp.
  • Gates where the cross correlation ratio is below min_rhi. Using the default parameter this filtering is disabled.
  • Gates where any of the above three fields are masked or contain invalid values (NaNs or infs).
  • If any of these three fields do not exist in the radar that fields filter criteria is not applied.
Parameters:
radar : Radar

Radar object from which the gate filter will be built.

refl_field, ncp_field, rhv_field : str

Names of the radar fields which contain the reflectivity, normalized coherent power (signal quality index) and cross correlation ratio (RhoHV) from which the gate filter will be created using the above criteria. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

min_ncp, min_rhv : float

Minimum values for the normalized coherence power and cross correlation ratio. Gates in these fields below these limits as well as gates which are masked or contain invalid values will be excluded and not used in calculation which use the filter. A value of None will disable filtering based upon the given field including removing masked or gates with an invalid value. To disable the thresholding but retain the masked and invalid filter set the parameter to a value below the lowest value in the field.

min_refl, max_refl : float

Minimum and maximum values for the reflectivity. Gates outside of this interval as well as gates which are masked or contain invalid values will be excluded and not used in calculation which use this filter. A value or None for one of these parameters will disable the minimum or maximum filtering but retain the other. A value of None for both of these values will disable all filtering based upon the reflectivity including removing masked or gates with an invalid value. To disable the interval filtering but retain the masked and invalid filter set the parameters to values above and below the lowest and greatest values in the reflectivity field.

Returns:
gatefilter : GateFilter

A gate filter based upon the described criteria. This can be used as a gatefilter parameter to various functions in pyart.correct.

pyart.correct.phase_proc_lp(radar, offset, debug=False, self_const=60000.0, low_z=10.0, high_z=53.0, min_phidp=0.01, min_ncp=0.5, min_rhv=0.8, fzl=4000.0, sys_phase=0.0, overide_sys_phase=False, nowrap=None, really_verbose=False, LP_solver='cylp', refl_field=None, ncp_field=None, rhv_field=None, phidp_field=None, kdp_field=None, unf_field=None, window_len=35, proc=1, coef=0.914)[source]

Phase process using a LP method [1].

Parameters:
radar : Radar

Input radar.

offset : float

Reflectivity offset in dBz.

debug : bool, optional

True to print debugging information.

self_const : float, optional

Self consistency factor.

low_z : float, optional

Low limit for reflectivity. Reflectivity below this value is set to this limit.

high_z : float, optional

High limit for reflectivity. Reflectivity above this value is set to this limit.

min_phidp : float, optional

Minimum Phi differential phase.

min_ncp : float, optional

Minimum normal coherent power.

min_rhv : float, optional

Minimum copolar coefficient.

fzl : float, optional

Maximum altitude.

sys_phase : float, optional

System phase in degrees.

overide_sys_phase : bool, optional

True to use sys_phase as the system phase. False will calculate a value automatically.

nowrap : int or None, optional

Gate number to begin phase unwrapping. None will unwrap all phases.

really_verbose : bool, optional

True to print LPX messaging. False to suppress.

LP_solver : ‘pyglpk’ or ‘cvxopt’, ‘cylp’, or ‘cylp_mp’, optional

Module to use to solve LP problem. Default is ‘cylp’.

refl_field, ncp_field, rhv_field, phidp_field, kdp_field : str, optional

Name of field in radar which contains the horizonal reflectivity, normal coherent power, copolar coefficient, differential phase shift, and differential phase. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

unf_field : str, optional

Name of field which will be added to the radar object which will contain the unfolded differential phase. Metadata for this field will be taken from the phidp_field. A value of None will use the default field name as defined in the Py-ART configuration file.

window_len : int, optional

Length of Sobel window applied to PhiDP field when prior to calculating KDP.

proc : int, optional

Number of worker processes, only used when LP_solver is ‘cylp_mp’.

coef : float, optional

Exponent linking Z to KDP in self consistency. kdp=(10**(0.1z))*coef

Returns:
reproc_phase : dict

Field dictionary containing processed differential phase shifts.

sob_kdp : dict

Field dictionary containing recalculated differential phases.

References

[1] Giangrande, S.E., R. McGraw, and L. Lei. An Application of Linear Programming to Polarimetric Radar Differential Phase Processing. J. Atmos. and Oceanic Tech, 2013, 30, 1716.

pyart.correct.phase_proc_lp_gf(radar, gatefilter=None, debug=False, self_const=60000.0, low_z=10.0, high_z=53.0, min_phidp=0.01, fzl=4000.0, system_phase=None, nowrap=None, really_verbose=False, LP_solver='cylp', refl_field=None, phidp_field=None, kdp_field=None, unf_field=None, window_len=35, proc=1, coef=0.914, ncpts=None, first_gate_sysp=None, offset=0.0, doc=0)[source]

Phase process using a LP method [1] using Py-ART’s Gatefilter.

Parameters:
radar : Radar

Input radar.

gatefilter : Gatefilter, optional

Py-ART gatefilter object indicating where processing should be carried out

debug : bool, optional

True to print debugging information.

self_const : float, optional

Self consistency factor.

low_z : float, optional

Low limit for reflectivity. Reflectivity below this value is set to this limit.

high_z : float, optional

High limit for reflectivity. Reflectivity above this value is set to this limit.

fzl : float, optional

Maximum altitude.

system_phase : float, optional

System phase in degrees.

nowrap : int or None, optional

Gate number to begin phase unwrapping. None will unwrap all phases.

really_verbose : bool, optional

True to print LPX messaging. False to suppress.

LP_solver : ‘pyglpk’ or ‘cvxopt’, ‘cylp’, or ‘cylp_mp’, optional

Module to use to solve LP problem. Default is ‘cylp’.

refl_field, ncp_field, rhv_field, phidp_field, kdp_field : str, optional

Name of field in radar which contains the horizonal reflectivity, normal coherent power, copolar coefficient, differential phase shift, and differential phase. A value of None for any of these parameters will use the default field name as defined in the Py-ART configuration file.

unf_field : str, optional

Name of field which will be added to the radar object which will contain the unfolded differential phase. Metadata for this field will be taken from the phidp_field. A value of None will use the default field name as defined in the Py-ART configuration file.

window_len : int, optional

Length of Sobel window applied to PhiDP field when prior to calculating KDP.

proc : int, optional

Number of worker processes, only used when LP_solver is ‘cylp_mp’.

coef : float, optional

Exponent linking Z to KDP in self consistency. kdp=(10**(0.1z))*coef

ncpts : int, optional

Minimum number of points in a ray. Regions within a ray smaller than this or beginning before this gate number are excluded from unfolding.

offset : float, optional

Reflectivity offset to add in dBz.

doc : int, optional

Number of gates to “doc” off the end of a ray.

Returns:
reproc_phase : dict

Field dictionary containing processed differential phase shifts.

sob_kdp : dict

Field dictionary containing recalculated differential phases.

References

[1] Giangrande, S.E., R. McGraw, and L. Lei. An Application of Linear Programming to Polarimetric Radar Differential Phase Processing. J. Atmos. and Oceanic Tech, 2013, 30, 1716.

pyart.correct.ptoa_to_sf(ptoa, pulse_width, wavelen, antenna_gain, coeff_band=1.2)[source]

Converts the sun power at the top of the atmosphere (in dBm) into solar flux.

Parameters:
ptoa : float

sun power at the top of the amosphere. It already takes into account the correction for antenna polarization

pulse_width : float

pulse width [s]

wavelen : float

radar wavelength [m]

antenna_gain : float

the antenna gain [dB]

coeff_band : float

multiplicative coefficient applied to the inverse of the pulse width to get the effective bandwidth

Returns:
s0 : float

solar flux [10e-22 W/(m2 Hz)]

References

Altube P., J. Bech, O. Argemi, T. Rigo, 2015: Quality Control of Antenna Alignment and Receiver Calibration Using the Sun: Adaptation to Midrange Weather Radar Observations at Low Elevation Angles

pyart.correct.scanning_losses(angle_step, beamwidth)[source]

Given the antenna beam width and the integration angle, compute the losses due to the fact that the sun is not a point target and the antenna is scanning

Parameters:
angle_step : float

integration angle [deg]

beamwidth : float

3 dB-beamwidth [deg]

Returns:
la : float

The losses due to the scanning of the antenna [dB positive]

References

Altube P., J. Bech, O. Argemi, T. Rigo, 2015: Quality Control of Antenna Alignment and Receiver Calibration Using the Sun: Adaptation to Midrange Weather Radar Observations at Low Elevation Angles

pyart.correct.selfconsistency_bias(radar, zdr_kdpzh_dict, min_rhohv=0.92, filter_rain=True, max_phidp=20.0, smooth_wind_len=5, doc=None, fzl=None, thickness=700.0, min_rcons=20, dphidp_min=2, dphidp_max=16, parametrization='None', refl_field=None, phidp_field=None, zdr_field=None, temp_field=None, iso0_field=None, hydro_field=None, rhohv_field=None, temp_ref='temperature', check_wet_radome=True, wet_radome_refl=25.0, wet_radome_len_min=4, wet_radome_len_max=8, wet_radome_ngates_min=180, valid_gates_only=False, keep_points=False, kdp_wind_len=12)[source]

Estimates reflectivity bias at each ray using the self-consistency algorithm by Gourley

Parameters:
radar : Radar

radar object

zdr_kdpzh_dict : dict

dictionary containing a look up table relating ZDR with KDP/Zh for different elevations

min_rhohv : float

minimum RhoHV value to consider the data valid

filter_rain : bool

If True the hydrometeor classification is going to be used to filter out all gates that are not rain

max_phidp : float

maximum PhiDP value to consider the data valid

smooth_wind_len : int

length of the smoothing window

doc : float

Number of gates at the end of each ray to to remove from the calculation.

fzl : float

Freezing layer, gates above this point are not included in the correction.

thickness : float

assumed melting layer thickness [m]

min_rcons : int

minimum number of consecutive gates to consider a valid segment of PhiDP

dphidp_min : float

minimum differential phase shift in a segment

dphidp_max : float

maximum differential phase shift in a segment

parametrization : str

The type of parametrization for the self-consistency curves. Can be ‘None’, ‘Gourley’, ‘Wolfensberger’, ‘Louf’, ‘Gorgucci’ or ‘Vaccarono’. ‘None’ will use tables contained in zdr_kdpzh_dict.

refl_field, phidp_field, zdr_field : str

Field names within the radar object which represent the reflectivity, differential phase and differential reflectivity fields. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_field, iso0_field, hydro_field, rhohv_field : str

Field names within the radar object which represent the temperature, the height relative to the iso0, the hydrometeor classification and the co-polar correlation fields. A value of None will use the default field name as defined in the Py-ART configuration file. They are going to be used only if available.

kdpsim_field, phidpsim_field : str

Field names which represent the estimated specific differential phase and differential phase. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_ref : str

the field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl

check_wet_radome : Bool

if True the average reflectivity of the closest gates to the radar is going to be check to find out whether there is rain over the radome. If there is rain no bias will be computed

wet_radome_refl : Float

Average reflectivity of the gates close to the radar to consider the radome as wet

wet_radome_len_min, wet_radome_len_max : int

Mim and max gate indices of the disk around the radome used to decide whether the radome is wet

wet_radome_ngates_min : int

Minimum number of valid gates to consider that the radome is wet

valid_gates_only : Bool

If True the reflectivity bias obtained for each valid ray is going to be assigned only to gates of the segment used. That will give more weight to longer segments when computing the total bias.

keep_points : Bool

If True the ZDR, ZH and KDP of the gates used in the self- consistency algorithm are going to be stored for further analysis

kdp_wind_len : int

The length of the window used to compute KDP with the single window least square method

Returns:
refl_bias_dict : dict

the bias at each ray field and metadata

selfconsistency_dict : dict

If keep_poinst set, a dictionary containing the measured valid values of ZDR, Zh and KDP. None otherwise

pyart.correct.selfconsistency_bias2(radar, zdr_kdpzh_dict, min_rhohv=0.92, min_zdr=0.2, filter_rain=True, max_phidp=20.0, smooth_wind_len=5, doc=None, fzl=None, thickness=700.0, min_rcons=20, parametrization='None', refl_field=None, phidp_field=None, zdr_field=None, temp_field=None, iso0_field=None, hydro_field=None, rhohv_field=None, kdp_field=None, temp_ref='temperature', check_wet_radome=True, wet_radome_refl=25.0, wet_radome_len_min=4, wet_radome_len_max=8, wet_radome_ngates_min=180, keep_points=False, bias_per_gate=False)[source]

Estimates reflectivity bias at each ray using the self-consistency algorithm by Gourley

Parameters:
radar : Radar

radar object

zdr_kdpzh_dict : dict

dictionary containing a look up table relating ZDR with KDP/Zh for different elevations

min_rhohv : float

minimum RhoHV value to consider the data valid

min_zdr : float

minimum ZDR value to consider the data valid

filter_rain : bool

If True the hydrometeor classification is going to be used to filter out all gates that are not rain

max_phidp : float

maximum PhiDP value to consider the data valid

smooth_wind_len : int

length of the smoothing window

doc : float

Number of gates at the end of each ray to to remove from the calculation.

fzl : float

Freezing layer, gates above this point are not included in the correction.

thickness : float

assumed melting layer thickness [m]

min_rcons : int

minimum number of consecutive gates to consider a valid segment of PhiDP

parametrization : str

The type of parametrization for the self-consistency curves. Can be ‘None’, ‘Gourley’, ‘Wolfensberger’, ‘Louf’, ‘Gorgucci’ or ‘Vaccarono’. ‘None’ will use tables contained in zdr_kdpzh_dict.

refl_field, kdp_field, zdr_field : str

Field names within the radar object which represent the reflectivity, differential phase and differential reflectivity fields. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_field, iso0_field, hydro_field, rhohv_field, phidp_field : str

Field names within the radar object which represent the temperature, the height relative to the iso0, the hydrometeor classification and the co-polar correlation fields. A value of None will use the default field name as defined in the Py-ART configuration file. They are going to be used only if available.

kdpsim_field, phidpsim_field : str

Field names which represent the estimated specific differential phase and differential phase. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_ref : str

the field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl

check_wet_radome : Bool

if True the average reflectivity of the closest gates to the radar is going to be check to find out whether there is rain over the radome. If there is rain no bias will be computed

wet_radome_refl : Float

Average reflectivity of the gates close to the radar to consider the radome as wet

wet_radome_len_min, wet_radome_len_max : int

Mim and max gate indices of the disk around the radome used to decide whether the radome is wet

wet_radome_ngates_min : int

Minimum number of valid gates to consider that the radome is wet

keep_points : Bool

If True the ZDR, ZH and KDP of the gates used in the self- consistency algorithm are going to be stored for further analysis

bias_per_gate : Bool

If True the bias per gate will be computed

Returns:
kdp_data_dict : dict

A dictionary containing valid observed and estimated using self- consistency values of KDP

refl_bias_dict : dict

If bias_per_gate is set, the bias at each gate field and metadata. None otherwise

selfconsistency_dict : dict

If keep_poinst set, a dictionary containing the measured valid values of ZDR, Zh and KDP. None otherwise

pyart.correct.selfconsistency_kdp_phidp(radar, zdr_kdpzh_dict, min_rhohv=0.92, filter_rain=True, max_phidp=20.0, smooth_wind_len=5, doc=None, fzl=None, thickness=700.0, parametrization='None', refl_field=None, phidp_field=None, zdr_field=None, temp_field=None, iso0_field=None, hydro_field=None, rhohv_field=None, kdpsim_field=None, phidpsim_field=None, temp_ref='temperature')[source]

Estimates KDP and PhiDP in rain from Zh and ZDR using a selfconsistency relation between ZDR, Zh and KDP. Private method

Parameters:
radar : Radar

radar object

zdr_kdpzh_dict : dict

dictionary containing a look up table relating ZDR with KDP/Zh for different elevations

min_rhohv : float

minimum RhoHV value to consider the data valid

filter_rain : bool

If True the hydrometeor classification is going to be used to filter out all gates that are not rain

max_phidp : float

maximum PhiDP value to consider the data valid

smooth_wind_len : int

length of the smoothing window

doc : float

Number of gates at the end of each ray to to remove from the calculation.

fzl : float

Freezing layer, gates above this point are not included in the correction.

thickness : float

assumed melting layer thickness [m]

parametrization : str

The type of parametrization for the self-consistency curves. Can be ‘None’, ‘Gourley’, ‘Wolfensberger’, ‘Louf’, ‘Gorgucci’ or ‘Vaccarono’. ‘None’ will use tables contained in zdr_kdpzh_dict.

refl_field, phidp_field, zdr_field : str

Field names within the radar object which represent the reflectivity, differential phase and differential reflectivity fields. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_field, iso0_field, hydro_field, rhohv_field : str

Field names within the radar object which represent the temperature, the height relative to the iso0, the hydrometeor classification and the co-polar correlation fields. A value of None will use the default field name as defined in the Py-ART configuration file. They are going to be used only if available.

kdpsim_field, phidpsim_field : str

Field names which represent the estimated specific differential phase and differential phase. A value of None will use the default field name as defined in the Py-ART configuration file.

temp_ref : str

the field use as reference for temperature. Can be either temperature, height_over_iso0 or fixed_fzl

Returns:
kdp_sim_dict, phidp_sim_dict : dict

the KDP and PhiDP estimated fields and metadata

pyart.correct.smooth_masked(raw_data, wind_len=11, min_valid=6, wind_type='median')[source]

smoothes the data using a rolling window. data with less than n valid points is masked.

Parameters:
raw_data : float masked array

The data to smooth.

window_len : float

Length of the moving window

min_valid : float

Minimum number of valid points for the smoothing to be valid

wind_type : str

type of window. Can be median or mean

Returns:
data_smooth : float masked array

smoothed data

pyart.correct.smooth_phidp_double_window(radar, ind_rmin=10, ind_rmax=500, min_rcons=11, zmin=20.0, zmax=40, swind_len=11, smin_valid=6, lwind_len=31, lmin_valid=16, zthr=40.0, psidp_field=None, refl_field=None, phidp_field=None)[source]

correction of the system offset and smoothing using two window

Parameters:
radar : Radar

Radar object for which to determine the system phase.

ind_rmin, ind_rmax : int

Min and max range index where to look for continuous precipitation

min_rcons : int

The minimum number of consecutive gates to consider it a rain cell.

zmin, zmax : float

Minimum and maximum reflectivity to consider it a rain cell

swind_len : int

Length of the short moving window used to smooth

smin_valid : int

Minimum number of valid bins to consider the short window smooth data valid

lwind_len : int

Length of the long moving window used to smooth

lmin_valid : int

Minimum number of valid bins to consider the long window smooth data valid

zthr : float

reflectivity value above which the short window is used

psidp_field : str

Field name within the radar object which represent the differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

refl_field : str

Field name within the radar object which represent the reflectivity. A value of None will use the default field name as defined in the Py-ART configuration file.

phidp_field : str

Field name within the radar object which represent the corrected differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

Returns:
phidp_dict : dict

The corrected phidp field

pyart.correct.smooth_phidp_single_window(radar, ind_rmin=10, ind_rmax=500, min_rcons=11, zmin=20.0, zmax=40, wind_len=11, min_valid=6, psidp_field=None, refl_field=None, phidp_field=None)[source]

correction of the system offset and smoothing using one window

Parameters:
radar : Radar

Radar object for which to determine the system phase.

ind_rmin, ind_rmax : int

Min and max range index where to look for continuous precipitation

min_rcons : int

The minimum number of consecutive gates to consider it a rain cell.

zmin, zmax : float

Minimum and maximum reflectivity to consider it a rain cell

wind_len : int

Length of the moving window used to smooth

min_valid : int

Minimum number of valid bins to consider the smooth data valid

psidp_field : str

Field name within the radar object which represent the differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

refl_field : str

Field name within the radar object which represent the reflectivity. A value of None will use the default field name as defined in the Py-ART configuration file.

phidp_field : str

Field name within the radar object which represent the corrected differential phase shift. A value of None will use the default field name as defined in the Py-ART configuration file.

Returns:
phidp_dict : dict

The corrected phidp field

pyart.correct.solar_flux_lookup(solar_flux, wavelen)[source]

Given the observed solar flux at 10.7 cm wavelength, returns the solar flux at the given radar wavelength

Parameters:
solar_flux : float array

the solar fluxes measured at 10.7 cm wavelength [10e-22 W/(m2 Hz)]

wavelen : float

radar wavelength [m]

Returns:
s0 : float

the radar flux at the radar wavelength [10e-22 W/(m2 Hz)]

References

Altube P., J. Bech, O. Argemi, T. Rigo, 2015: Quality Control of Antenna Alignment and Receiver Calibration Using the Sun: Adaptation to Midrange Weather Radar Observations at Low Elevation Angles

pyart.correct.sun_position_mfr(dt, lat_deg, lon_deg, refraction=True)[source]

Calculate the sun position for the given time (dt) at the given position (lat, lon).

Parameters:
dt : datetime object

the time when to look for the sun

lat_deg, lon_deg: floats

latitude and longitude in degrees

refraction : boolean

whether to correct for refraction or not

Returns:
elev_sun, azim_sun : floats

elevation and azimuth angles of the sun respect to the sensor in degrees

pyart.correct.sun_position_pysolar(dt, lat, lon, elevation=0.0)[source]

obtains the sun position in antenna coordinates using the pysolar library.

Parameters:
dt : datetime object

the time when to look for the sun

lat, lon : float

latitude and longitude of the sensor in degrees

Returns:
el, az : float

elevation and azimuth angles of the sun respect to the sensor in degrees

pyart.correct.sun_power(solar_flux, pulse_width, wavelen, antenna_gain, angle_step, beamwidth, coeff_band=1.2)[source]

computes the theoretical sun power detected at the antenna [dBm] as it would be without atmospheric attenuation (sun power at top of the atmosphere) for a given solar flux and radar characteristics

Parameters:
solar_flux : float array

the solar fluxes measured at 10.7 cm wavelength [10e-22 W/(m2 Hz)]

pulse_width : float

pulse width [s]

wavelen : float

radar wavelength [m]

antenna_gain : float

the antenna gain [dB]

angle_step : float

integration angle [deg]

beamwidth : float

3 dB-beamwidth [deg]

coeff_band : float

multiplicative coefficient applied to the inverse of the pulse width to get the effective bandwidth

Returns:
pwr_det : float array

the detected power

References

Altube P., J. Bech, O. Argemi, T. Rigo, 2015: Quality Control of Antenna Alignment and Receiver Calibration Using the Sun: Adaptation to Midrange Weather Radar Observations at Low Elevation Angles

pyart.correct.sun_retrieval(az_rad, az_sun, el_rad, el_sun, sun_hit, sun_hit_std, az_width_co=None, el_width_co=None, az_width_cross=None, el_width_cross=None, is_zdr=False)[source]

Estimates sun parameters from sun hits

Parameters:
az_rad, az_sun, el_rad, el_sun : float array

azimuth and elevation values of the sun and the radar

sun_hit : float array

sun hit value. Either power in dBm or ZDR in dB

sun_hit_std : float array

standard deviation of the sun hit value in dB

az_width_co, el_width_co, az_width_cross, el_width_cross : float

azimuth and elevation antenna width for each channel

is_zdr : boolean

boolean to signal that is ZDR data

Returns:
val, val_std : float

retrieved value and its standard deviation

az_bias, el_bias : float

retrieved azimuth and elevation antenna bias respect to the sun position

az_width, el_width : float

retrieved azimuth and elevation antenna widths

nhits : int

number of sun hits used in the retrieval

par : float array

and array with the 5 parameters of the Gaussian fit