SpectraProcessing

This routine calculates the radar moments for the RPG 94 GHz FMCW radar ‘LIMRAD94’ and generates a NetCDF4 file. The generated files can be used as input for the Cloudnet processing chain.

param **date:format YYYYMMDD
type **date:string
param **path:path where NetCDF file will be stored
type **path:string

Example

python spec2mom_limrad94.py date=20181201 path=/tmp/pycharm_project_626/scripts_Willi/cloudnet_input/
pyLARDA.SpectraProcessing.replace_fill_value(data, newfill)[source]

Replaces the fill value of an spectrum container by their time and range specific mean noise level.

Parameters:
  • data (numpy.array) – 3D spectrum array (time, range, velocity)
  • newfill (numpy.array) – 2D new fill values for 3rd dimension (velocity)
Returns:

spectrum with mean noise

Return type:

var (numpy.array)

pyLARDA.SpectraProcessing.estimate_noise_hs74[source]

//github.com/ARM-DOE/pyart/blob/master/pyart/util/hildebrand_sekhon.py

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

Parameters:
  • spectrum (array) – Doppler spectrum in linear units.
  • navg (int, optional) – The number of spectral bins over which a moving average has been taken. Corresponds to the p variable from equation 9 of the article. The default value of 1 is appropriate when no moving average has been applied to the spectrum.
  • std_div (float, optional) – Number of standard deviations above mean noise floor to specify the signal threshold, default: threshold=mean_noise + 6*std(mean_noise)
  • nnoise_min (int, optional) – Minimum number of noise samples to consider the estimation valid.
Returns:

tuple with

  • mean (float): Mean of points in the spectrum identified as noise.
  • threshold (float): Threshold separating noise from signal. The point in the spectrum with this value or below should be considered as noise, above this value signal. It is possible that all points in the spectrum are identified as noise. If a peak is required for moment calculation then the point with this value should be considered as signal.
  • var (float): Variance of the points in the spectrum identified as noise.
  • nnoise (int): Number of noise points in the spectrum.

References

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

Type:REFERENCE TO ARM PYART GITHUB REPO
Type:https
pyLARDA.SpectraProcessing.find_peak_edges[source]

Returns the indices of left and right edge of the main signal peak in a Doppler spectra.

Parameters:
  • signal (numpy.array) – 1D array Doppler spectra
  • threshold – noise threshold
Returns:

indices of signal minimum/maximum velocity

Return type:

[index_left, index_right] (list)

pyLARDA.SpectraProcessing.radar_moment_calculation[source]

reflectivity, mean Doppler velocity, spectral width, skewness, and kurtosis of one Doppler spectrum. Optimized for the use of Numba.

Note

Divide the signal_sum by 2 because vertical and horizontal channel are added. Subtract half of of the Doppler resolution from mean Doppler velocity, because

Parameters:
  • signal (-) – detected signal from a Doppler spectrum
  • vel_bins (-) – extracted velocity bins of the signal (same length as signal)
  • DoppRes (-) – resolution of the Doppler spectra (different for each chirp)
Returns:

dict containing

  • Ze_lin (float array): reflectivity (0.Mom) over range of velocity bins [mm6/m3]
  • VEL (float array): mean velocity (1.Mom) over range of velocity bins [m/s]
  • sw (float array): spectrum width (2.Mom) over range of velocity bins [m/s]
  • skew (float array): skewness (3.Mom) over range of velocity bins
  • kurt (float array): kurtosis (4.Mom) over range of velocity bins

Type:

Calculation of radar moments

pyLARDA.SpectraProcessing.despeckle[source]

Remove small patches (speckle) from any given mask by checking 5x5 box around each pixel, more than half of the points in the box need to be 1 to keep the 1 at current pixel

Parameters:
  • mask (numpy.array, integer) – 2D mask where 1 = an invalid/fill value and 0 = a data point (time, height)
  • min_percentage (float) – minimum percentage of neighbours that need to be signal above noise
Returns:

mask … speckle-filtered matrix of 0 and 1 that represents (cloud) mask [height x time]

pyLARDA.SpectraProcessing.make_container_from_spectra(spectra_all_chirps, values, paraminfo, invalid_mask, varname='')[source]

This routine will generate a larda container from calculated moments from spectra.

Parameters:
  • spectra_all_chirps (list of dicts) – dimension [nchirps], containing the spectrum values of the 94 GHz RPG cloud radar
  • values (numpy array) – dimension [nrange, ntimes], values of calculated moments
  • paraminfo (dict) – information from params_[campaign].toml for the specific variable
Returns:

larda data container

Return type:

container (dict)

pyLARDA.SpectraProcessing.load_spectra_rpgfmcw94(larda, time_span, rpg_radar='LIMRAD94', **kwargs)[source]

This routine will generate a list of larda containers including spectra of the RPG-FMCW 94GHz radar. The list-container at return will contain the additional information, for each chirp.

Parameters:
  • rpg_radar (string) – name of the radar system as defined in the toml file
  • larda (class larda) – Initialized pyLARDA, already connected to a specific campaign
  • time_span (list) – Starting and ending time point in datetime format.
  • **noise_factor (float) – Noise factor, number of standard deviations from mean noise floor
  • **ghost_echo_1 (bool) – Filters ghost echos which occur over all chirps during precipitation.
  • **ghost_echo_2 (bool) – Filters ghost echos which occur over 1 chirps during precipitation.
  • **estimate_noise (boal) –

    If True, adds the following noise estimation values to the container:

    • mean (2d ndarray): Mean noise level of the spectra.
    • threshold (2d ndarray): Noise threshold, values above this threshold are consider as signal.
    • variance (2d ndarray): The variance of the mean noise level.
    • numnoise (2d ndarray): Number of Pixels that are cconsideras noise.
    • signal (2d ndarray): Boolean array, a value is True if no signal was detected.
    • bounds (3d ndarrax): Dimensions [n_time, n_range, 2] containing the integration boundaries.
Returns:

list of larda data container

  • spec[i_chirps][‘no_av’] (float): Number of spectral averages divided by the number of FFT points
  • spec[i_chirps][‘DoppRes’] (float): Doppler resolution for
  • spec[i_chirps][‘SL’] (2D-float): Sensitivity limit (dimensions: time, range)
  • spec[i_chirps][‘NF’] (string): Noise factor, default = 6.0
  • spec[i_chirps][‘rg_offsets’] (list): Indices, where chipr shifts

Return type:

container (list)

pyLARDA.SpectraProcessing.dealiasing_check(masked3D)[source]

Checks for folding.

Parameters:
  • masked3D (numpy.array) – 3D (time, range, velocity)
  • vel (list) – contains 1D numpy.arrays for each chirp
  • mean_noise (numpy.array) –
  • rg_offsets (numpy.array) –
Returns:

Return type:

alias_flag (numpy.array)

pyLARDA.SpectraProcessing.dealiasing(spectra: numpy.array, vel_bins_per_chirp: List[numpy.array], noisefloor: numpy.array, rg_offsets: List[int] = None, show_triple: bool = False, vel_offsets: List[int] = None, jump: int = None) → Union[numpy.array, List[numpy.array]][source]

Peaks exceeding the maximum unambiguous Doppler velocity range of ± v_Nyq in [m s-1] appear at the next upper (lower) range gate at the other end of the velocity spectrum. The dealiasing method presented here aims to correct for this and is applied to every time step.

Logging level INFO shows the de-aliasing progress bar.

Parameters:
  • spectra – dim = (n_time, n_range, n_velocity) in linear units!
  • vel_bins_per_chirp – len = (n_chirp), each list element contains a numpy array of velocity bins
  • noisefloor – dim = (n_time, n_range) in linear units!
  • rg_offsets (optional) – dim = (n_chirp + 1), starting with 0, range indices where chirp shift occurs
  • show_triple (optional) – if True, return dealiased spectra including the triplication
  • vel_offsets (optional) – velocity window around the main peak [x1, x2], x1 < 0, x2 > 0 ! in [m s-1], default [-6.0, +9.0]
  • jump (optional) – maximum number of Doppler bins a spectrum can change in two adjacent range bins
Returns:

tuple containing

  • dealiased_spectra: dim = (n_time, n_range, 3 * n_velocity), de-aliased Doppler spectrum
  • dealiased_mask: dim = (n_time, n_range, 3 * n_velocity), True if no signal
  • velocity_new: len = (n_chirp), each list element contains a numpy array of velocity bins for the respective chirp of ± 3*v_Nyq in [m s-1]
  • signal_boundaries: indices of left and right edge of a signal, [-1, -1] if no signal
  • search_path: indices of left and right edge of the search path, [-1, -1] if no signal
  • idx_peak_matrix: indices of the main peaks, [NDbins / 2] if no signal

Todo

  • add time–height mask for dealiasing
  • search window relative to [m s-1]
  • abs(idx_new_peak - mean_idx_last_ts) > 120: –> relativ
pyLARDA.SpectraProcessing.noise_estimation_uncompressed_data(data, n_std=6.0, **kwargs)[source]

Creates a dict containing the noise threshold, mean noise level, the variance of the noise, the number of noise values in the spectrum, and the boundaries of the main signal peak, if there is one

Parameters:
  • data (dict) – data container, containing data[‘var’] of dimension (n_ts, n_range, n_Doppler_bins)
  • **n_std_deviations (float) – threshold = number of standard deviations above mean noise floor, default: threshold is the value of the first non-noise value
Returns:

dict with noise floor estimation for all time and range points

pyLARDA.SpectraProcessing.mira_noise_calculation(radar_const, SNRCorFaCo, HSDco, noise_power_co, range_ka)[source]
Parameters:
  • radar_const
  • SNRCorFaCo
  • HSDco
  • noise_power_co
  • range_ka
Returns:

noise level in linear units

pyLARDA.SpectraProcessing.despeckle2D(data, min_perc=80.0)[source]

This function is used to remove all spectral lines for one time-range-pixel if surrounding% of the sourounding pixels are fill_values.

Parameters:data (numpy.array) – cloud radar Doppler spectra, dimensions: (time, range, velocity), unit: [mm6 mm-3 m s-1]
Keyword Arguments:
 min_perc (float) – minimum percentage value of neighbouring pixel, that need to be above the noise threshold
Returns:where True = fill_value, and False = signal, dimensions: (time, range, velocity)
Return type:mask (numpy.array, bool)
pyLARDA.SpectraProcessing.filter_ghost_1(data, rg, vel, offset, dBZ_thresh=-20.0, reduce_by=1.5, **kwargs)[source]

This function is used to remove certain spectral lines “speckle ghost echoes” from all chirps of RPG FMCW 94GHz cloud radar spectra. The speckle occur usually near the maximum unambiguous Doppler velocity.

Parameters:
  • data (numpy.array) – cloud radar Doppler spectra, dimensions: (time, range, velocity), unit: [mm6 mm-3 m s-1]
  • rg (numpy.array) – range values, unit [m]
  • vel (list of numpy.arrays) – contains the Doppler velocity values for each chirp, dimension = n_chirps
  • offset (list, integer) – range indices where the chirp changes takes place, dimension = n_chirps + 1 (starting with 0)
  • dBZ_thresh (float) – values below will be considered as ghost echo
  • reduce_by (float) – reduce the maximum unambiguous Doppler velocity by this amount in [m s-1]
  • **ignore_chirp1 (bool) – Don’t filter ghost echos of this type for first chirp (set to True if not given)
  • **Z_thresh (float) – Ze in dBZ to be exceeded in lowest 500 m range for filter to be activated
Returns:

where True = fill_value, and False = signal, dimensions: (time, range, velocity)

Return type:

mask (numpy.array, bool)

pyLARDA.SpectraProcessing.filter_ghost_2(data, rg, SL, first_offset, dBZ_thresh=-5.0, reduce_by=10.0)[source]

This function is used to remove curtain-like ghost echoes from the first chirp of RPG FMCW 94GHz cloud radar spectra.

Parameters:
  • data (numpy.array) – cloud radar Doppler spectra, dimensions: (time, range, velocity), unit: [mm6 mm-3 m s-1]
  • rg (numpy.array) – range values, unit [m]
  • SL (numpy.array) – sensitivity limit, dimension: (time, range), unit: [mm6 mm-3]
  • first_offset (integer) – range index where the first chirp change takes place
  • dBZ_thresh (float) – minimum threshold in [dBZ], where ghost echos can be assumed
  • reduce_by (float) – reduce the sensitivity limit by this amount of [dBZ]
Returns:

where True = fill_value, and False = signal, dimensions: (time, range, velocity)

Return type:

mask (numpy.array, bool)

pyLARDA.SpectraProcessing.spectra2moments(ZSpec, paraminfo, **kwargs)[source]

This routine calculates the radar moments: reflectivity, mean Doppler velocity, spectrum width, skewness and kurtosis from the level 0 spectrum files of the 94 GHz RPG cloud radar.

Parameters:
  • ZSpec (dict) – list containing the dicts for each chrip of RPG-FMCW Doppler cloud radar
  • paraminfo (dict) – information from params_[campaign].toml for the system LIMRAD94
Returns:

dictionary of larda containers, including larda container for Ze, VEL, sw, skew, kurt

Return type:

container_dict (dict)

pyLARDA.SpectraProcessing.heave_correction(moments, date, path_to_seapath='/projekt2/remsens/data_new/site-campaign/rv_meteor-eurec4a/instruments/RV-METEOR_DSHIP', mean_hr=True, only_heave=False, use_cross_product=True, transform_to_earth=True, add=False)[source]

Correct mean Doppler velocity for heave motion of ship (RV-Meteor) Calculate heave rate from seapath measurements and create heave correction array. If Doppler velocity is given as an input, correct it and return an array with the corrected Doppler velocities. Without Doppler Velocity input, only the heave correction array is returned.

Parameters:
  • moments – LIMRAD94 moments container as returned by spectra2moments in spec2mom_limrad94.py, C1/2/3_Range, SeqIntTime and Inc_ElA (for time (ts)) from LV1 file
  • date (datetime.datetime) – object with date of current file
  • path_to_seapath (string) – path where seapath measurement files (daily dat files) are stored
  • mean_hr (bool) – whether to use the mean heave rate over the SeqIntTime or the heave rate at the start time of the chirp
  • only_heave (bool) – whether to use only heave to calculate the heave rate or include pitch and roll induced heave
  • use_cross_product (bool) – whether to use the cross product like Hannes Griesche https://doi.org/10.5194/amt-2019-434
  • transform_to_earth (bool) – transform cross product to earth coordinate system as described in https://repository.library.noaa.gov/view/noaa/17400
  • add (bool) – whether to add the heave rate or subtract it
Returns:

A number of variables

  • new_vel (ndarray); corrected Doppler velocities, same shape as moments[“VEL”][“var”] or list if no Doppler Velocity is given;
  • heave_corr (ndarray): heave rate closest to each radar timestep for each height bin, same shape as moments[“VEL”][“var”];
  • seapath_out (pd.DataFrame): data frame with all heave information from the closest time steps to the chirps

pyLARDA.SpectraProcessing.heave_correction_spectra(data, date, path_to_seapath='/projekt2/remsens/data_new/site-campaign/rv_meteor-eurec4a/instruments/RV-METEOR_DSHIP', mean_hr=True, only_heave=False, use_cross_product=True, transform_to_earth=True, add=False, **kwargs)[source]

Shift Doppler spectra to correct for heave motion of ship (RV-Meteor) Calculate heave rate from seapath measurements and create heave correction array. Translate the heave correction to a number spectra bins by which to move each spectra. If Spectra are given, shift them and return a 3D array with the shifted spectra. Without spectra input, only the heave correction array and the array with the number if bins to move is returned.

Parameters:
  • data – LIMRAD94 data container filled with spectra and C1/2/3_Range, SeqIntTime, MaxVel, DoppLen from LV1 file; for Claudia’s version the mean Doppler velocity is also needed
  • date (datetime.datetime) – object with date of current file
  • path_to_seapath (string) – path where seapath measurement files (daily dat files) are stored
  • mean_hr (bool) – whether to use the mean heave rate over the SeqIntTime or the heave rate at the start time of the chirp
  • only_heave (bool) – whether to use only heave to calculate the heave rate or include pitch and roll induced heave
  • use_cross_product (bool) – whether to use the cross product like Hannes Griesche https://doi.org/10.5194/amt-2019-434
  • transform_to_earth (bool) – transform cross product to earth coordinate system as described in https://repository.library.noaa.gov/view/noaa/17400
  • add (bool) – whether to add the heave rate or subtract it
  • **kwargs – shift (int): number of time steps to shift seapath data version (str): which version to use, ‘ca’ or ‘jr’
Returns:

A number of variables

  • new_spectra (ndarray); corrected Doppler velocities, same shape as data[“VHSpec”][“var”] or list if no Doppler Spectra are given;
  • heave_corr (ndarray): heave rate closest to each radar timestep for each height bin, shape = (time x range);
  • seapath_out (pd.DataFrame): data frame with all heave information from the closest time steps to the chirps

pyLARDA.SpectraProcessing.read_seapath(date, path='/projekt2/remsens/data_new/site-campaign/rv_meteor-eurec4a/instruments/RV-METEOR_DSHIP', **kwargs)[source]

Read in daily Seapath measurements from RV Meteor from .dat files to a pandas.DataFrame :param date: object with date of current file :type date: datetime.datetime :param path: path to seapath files :type path: str :param kwargs for read_csv: output_format (str): whether a pandas data frame or a xarray dataset is returned

Returns:DataFrame with Seapath measurements
Return type:seapath (DataFrame)
pyLARDA.SpectraProcessing.read_dship(date, **kwargs)[source]

Read in 1 Hz DSHIP data and return pandas DataFrame

Parameters:
Returns:

pd.DataFrame with 1 Hz DSHIP data

pyLARDA.SpectraProcessing.f_shiftTimeDataset(dataset)[source]

author: Claudia Acquistapace date: 25 november 2020 goal : function to shift time variable of the dataset to the central value of the time interval of the time step input:

dataset: xarray dataset
output:
dataset: xarray dataset with the time coordinate shifted added to the coordinates and the variables now referring to the shifted time array
pyLARDA.SpectraProcessing.calc_heave_rate(seapath, x_radar=-11, y_radar=4.07, z_radar=15.8, only_heave=False, use_cross_product=True, transform_to_earth=True)[source]

Calculate heave rate at a certain location of a ship with the measurements of the INS

Parameters:
  • seapath (pd.DataFrame) – Data frame with heading, roll, pitch and heave as columns
  • x_radar (float) – x position of location with respect to INS in meters
  • y_radar (float) – y position of location with respect to INS in meters
  • z_radar (float) – z position of location with respect to INS in meters
  • only_heave (bool) – whether to use only heave to calculate the heave rate or include pitch and roll induced heave
  • use_cross_product (bool) – whether to use the cross product like Hannes Griesche https://doi.org/10.5194/amt-2019-434
  • transform_to_earth (bool) – transform cross product to earth coordinate system as described in https://repository.library.noaa.gov/view/noaa/17400
Returns:

Data frame as input with additional columns radar_heave, pitch_heave, roll_heave and “heave_rate”

Return type:

seapath (pd.DataFrame)

pyLARDA.SpectraProcessing.f_calcRMatrix(rollShipArr, pitchShipArr, yawShipArr, NtimeShip)[source]

author: Claudia Acquistapace date : 27/10/2020 goal: function to calculate R matrix given roll, pitch, yaw input:

rollShipArr: roll array in degrees pitchShipArr: pitch array in degrees yawShipArr: yaw array in degrees NtimeShip: dimension of time array for the definition of R_inv as [3,3,dimTime]
output:
R[3,3,Dimtime]: array of rotational matrices, one for each time stamp
pyLARDA.SpectraProcessing.calc_heave_rate_claudia(data, x_radar=-11, y_radar=4.07, z_radar=-15.8)[source]

Calculate heave rate at a certain location on a ship according to Claudia Acquistapace’s approach

Parameters:
  • data (xr.DataSet) – Data Set with heading, roll, pitch and heave as columns
  • x_radar (float) – x position of location with respect to INS in meters
  • y_radar (float) – y position of location with respect to INS in meters
  • z_radar (float) – z position of location with respect to INS in meters

Returns: xr.DataSet with additional variable heave_rate

pyLARDA.SpectraProcessing.find_mdv_time_series(mdv_values, radar_time, n_ts_run)[source]

author: Claudia Acquistapace, Johannes Roettenbacher Identify, given a mean doppler velocity matrix, a sequence of length n_ts_run of values in the matrix at a given height that contains the minimum possible amount of nan values in it.

Parameters:
  • mdv_values (ndarray) – time x heigth matrix of Doppler Velocity
  • radar_time (ndarray) – corresponding radar time stamps in seconds (unix time)
  • n_ts_run (int) – number of timestamps needed in a mdv series
Returns:

time series of Doppler velocity with length n_ts_run time_series (ndarray): corresponding time stamps to Doppler velocity time series i_height_sel (int): index of chosen height valuesColumnMean (ndarray): time series of mean Doppler velocity averaged over height with length n_ts_run

Return type:

valuesTimeSerie (ndarray)

pyLARDA.SpectraProcessing.calc_time_shift(w_radar_meanCol, delta_t_min, delta_t_max, resolution, w_ship_chirp, timeSerieRadar, pathFig, chirp, hour, date)[source]

author: Claudia Acquistapace, Jan. H. Schween, Johannes Roettenbacher goal: calculate and estimation of the time lag between the radar time stamps and the ship time stamp

NOTE: adding or subtracting the obtained time shift depends on what you did during the calculation of the covariances: if you added/subtracted time _shift to t_radar you have to do the same for the ‘exact time’ Here is the time shift analysis as plot: <ww> is short for <w’_ship*w’_radar> i.e. covariance between vertical speeds from ship movements and radar its maximum gives an estimate for optimal agreement in vertical velocities of ship and radar <Delta w^2> is short for <(w[i]-w[i-1])^2> where w = w_rad - 2*w_ship - this is a measure for the stripeness. Its minimum gives an estimate how to get the smoothest w data :param w_radar_meanCol: time series of mean Doppler velocity averaged over height with no nan values :type w_radar_meanCol: ndarray :param delta_t_min: minimum time shift :type delta_t_min: float :param delta_t_max: maximum time shift :type delta_t_max: float :param resolution: time step by which to increment possible time shift :type resolution: float :param w_ship_chirp: vertical velocity of the radar at the exact chirp time step :type w_ship_chirp: ndarray :param timeSerieRadar: time stamps of the mean Doppler velocity time series (w_radar_meanCol) :type timeSerieRadar: ndarray :param pathFig: file path where figures should be stored :type pathFig: str :param chirp: which chirp is being processed :type chirp: int :param hour: which hour of the day is being processed (0-23) :type hour: int :param date: which day is being processed :type date: datetime

Returns: time shift between radar data and ship data in seconds, quicklooks for each calculation

pyLARDA.SpectraProcessing.calc_chirp_timestamps(radar_ts, date, version)[source]

Calculate the exact timestamp for each chirp corresponding with the center or start of the chirp The timestamp in the radar file corresponds to the end of a chirp sequence with an accuracy of 0.1 s

Parameters:
  • radar_ts (ndarray) – timestamps of the radar with milliseconds in seconds
  • date (datetime.datetime) – date which is being processed
  • version (str) – should the timestamp correspond to the ‘center’ or the ‘start’ of the chirp

Returns: dict with chirp timestamps

pyLARDA.SpectraProcessing.calc_shifted_chirp_timestamps(radar_ts, radar_mdv, chirp_ts, rg_borders_id, n_ts_run, Cs_w_radar, **kwargs)[source]

Calculates the time shift between each chirp time stamp and the ship time stamp for every hour and every chirp. :param radar_ts: radar time stamps in seconds (unix time) :type radar_ts: ndarray :param radar_mdv: time x height matrix of mean Doppler velocity from radar :type radar_mdv: ndarray :param chirp_ts: exact chirp time stamps :type chirp_ts: ndarray :param rg_borders_id: indices of chirp boundaries :type rg_borders_id: ndarray :param n_ts_run: number of time steps necessary for mean Doppler velocity time series :type n_ts_run: int :param Cs_w_radar: function of vertical velocity of radar against time :type Cs_w_radar: scipy.interpolate.CubicSpline :param **kwargs: no_chirps (int): number of chirps in radar measurement

plot_fig (bool): plot quicklook

Returns: time shifted chirp time stamps, array with time shifts for each chirp and hour

pyLARDA.SpectraProcessing.calc_corr_matrix_claudia(radar_ts, radar_rg, rg_borders_id, chirp_ts_shifted, Cs_w_radar)[source]

Calculate the correction matrix to correct the mean Doppler velocity for the ship vertical motion. :param radar_ts: original radar time stamps in seconds (unix time) :type radar_ts: ndarray :param radar_rg: radar range gates :type radar_rg: ndarray :param rg_borders_id: indices of chirp boundaries :type rg_borders_id: ndarray :param chirp_ts_shifted: hourly shifted chirp time stamps :type chirp_ts_shifted: dict :param Cs_w_radar: function of vertical velocity of radar against time :type Cs_w_radar: scipy.interpolate.CubicSpline

Returns: correction matrix for mean Doppler velocity

pyLARDA.SpectraProcessing.get_range_bin_borders(no_chirps, container)[source]

get the range bins which correspond to the chirp borders of a FMCW radar

Parameters:
  • no_chirps (int) – Number of chirps
  • container (dict) – Dictionary with C1/2/3Range variable from LV1 files
Returns:

ndarray with chirp borders including 0 range_bins

pyLARDA.SpectraProcessing.calc_heave_corr(container, chirp_ts, seapath, mean_hr=True)[source]

Calculate heave correction for mean Doppler velocity

Parameters:
  • container (larda container) – LIMRAD94 C1/2/3_Range, SeqIntTime, ts
  • chirp_ts (dict) – dictionary with exact radar chirp time stamps
  • seapath (pd.DataFrame) – Data frame with heave rate column (“heave_rate_radar”)
  • mean_hr (bool) – whether to use the mean heave rate over the SeqIntTime or the heave rate at the start time of the chirp
Returns:

heave rate closest to each radar timestep for each height bin, time x range

Return type:

heave_corr (ndarray)

pyLARDA.SpectraProcessing.calc_dopp_res(MaxVel, DoppLen, no_chirps, range_bins)[source]
Parameters:
  • MaxVel (ndarray) – Unambiguous Doppler velocity (+/-) m/s from LV1 file
  • DoppLen (ndarray) – Number of spectral lines in Doppler spectra from LV1 file
  • no_chirps (int) – Number of chirps
  • range_bins (ndarray) – range bin number of lower chirp borders, starts with 0
Returns:

1D array with Doppler resolution for each height bin

pyLARDA.SpectraProcessing.heave_rate_to_spectra_bins(heave_corr, doppler_res)[source]

translate the heave correction to Doppler spectra bins

Parameters:
  • heave_corr (ndarray) – heave rate closest to each radar timestep for each height bin, time x range
  • doppler_res (ndarray) – Doppler resolution of each chirp of LIMRAD94 for whole range 1 x range
Returns:

ndarray with number of bins to move each Doppler spectrum

  • n_dopp_bins_shift (ndarray): of same dimension as heave_corr
  • heave_corr

pyLARDA.SpectraProcessing.shift_seapath(seapath, shift)[source]

Shift seapath values by given shift

Parameters:
  • seapath (pd.Dataframe) – Dataframe with heave motion of RV-Meteor
  • shift (int) – number of time steps to shift data
Returns:

shifted Dataframe

pyLARDA.SpectraProcessing.find_closest_timesteps(df, ts)[source]

Find closest time steps in a dataframe to a time series

Parameters:
  • df (pd.DataFrame) – DataFrame with DatetimeIndex
  • ts (ndarray) – array with time stamps in unix format (seconds since 1-1-1970)
Returns:

pd.DataFrame with only the closest time steps to ts

pyLARDA.SpectraProcessing.spectra2sldr(ZSpec, paraminfo, **kwargs)[source]

This routine calculates the

Parameters:
  • ZSpec (dict) – list containing the dicts for each chrip of RPG-FMCW Doppler cloud radar
  • paraminfo (dict) – information from params_[campaign].toml for the system LIMRAD94
Returns:

dictionary of larda containers, including larda container for Ze, VEL, sw, skew, kurt

Return type:

container_dict (dict)

pyLARDA.SpectraProcessing.spectra2polarimetry(ZSpec, paraminfo, **kwargs)[source]

This routine calculates the

Parameters:
  • ZSpec (dict) – list containing the dicts for each chrip of RPG-FMCW Doppler cloud radar
  • paraminfo (dict) – information from params_[campaign].toml for the system LIMRAD94
Returns:

dictionary of larda containers, including larda container for Ze, VEL, sw, skew, kurt

Return type:

container_dict (dict)