lidalign package
Submodules
lidalign.SSC module
This module contains functions for the Sea-Surfac-Leveling (following Rott et al. 2022) and Sea-Surface-Calibration (following Meyer et al. 2026) from volumetric scans of the surrounding water surface. Many parts of this module are directly based on the published scripts by Rott et. al 2021 (see https://zenodo.org/record/5654866).
SSL: Used to align the scanning lidar offshore (reduce roll and pitch as much as possible) SSC: Calibration for scanning lidar elevation offset determination (determine offset and correct for afterwards in the software)
Authors: Paul Meyer, Andreas Rott Date: 2026-06-01
- class lidalign.SSC.CNRFitResult(success: bool, r_water: float, residuals: float = inf, x: list = None, params: dict = None, fit_data: Dataset = None, fit: Figure = None)[source]
Bases:
objectresult dataclass from CNR Sigmoid fits
- fit: Figure = None
- fit_data: Dataset = None
- params: dict = None
- r_water: float
- residuals: float = inf
- success: bool
- x: list = None
- class lidalign.SSC.EarthCurvature[source]
Bases:
object- static get_height(distance: float) float[source]
Get the height difference due to earth curvature for a given distance
Following Gramitzky, K., Jäger, F., Callies, D., Hildebrand, T., Lundquist, J. K., & Pauscher, L. (2025). Alignment of Scanning Lidars in Offshore Campaigns – an Extension of the Sea Surface Levelling Method. Wind and the atmosphere/Wind and turbulence. https://doi.org/10.5194/wes-2025-191
- Parameters:
distance (float) – distance from the observer [m]
- Returns:
height difference due to earth curvature [m]
- Return type:
float
- static get_intercept_with_curvature(height: float, elevation: float, max_distance: float = 10000) float[source]
Get interception point between a line (beam of lidar) and a circle (earth curvature) This is required if you want to find out, where your line-of-sight would hit the water.
- Parameters:
height (float) – height of the observer (lidar lense) in m
elevation (float) – elevation angle of the view (line-of-sight), elevation is defined positive up. Must be large enough, that the water can be seen at this elevation angle at all.
- Returns:
X-value where line of sight and circle have their first crossing (positive x only)
- Return type:
float
- class lidalign.SSC.GaussianTruncatedPulse(FWHM: float = 100, FWHM_Width_ratio: float = 2.6)[source]
Bases:
PulseShapePulse Shape of a truncated Gaussian Pulse. Child class of PulseShape.
Pulse shape definition following Schlipf, D. (2015). PhD Thesis: Lidar-Assisted Control Concepts for Wind Turbines.
- FWHM: float = 100
- FWHM_Width_ratio: float = 2.6
- static fit_weighting_to_data(r_data: ndarray[tuple[Any, ...], dtype[float64]], signal_strength: ndarray[tuple[Any, ...], dtype[float64]]) ndarray[tuple[Any, ...], dtype[float64]][source]
Fit the Pulse shape to given radial distribution of Signal strength. Can be used e.g., to fit the truncated gaussian to the CNR signal of a hard target.
- Parameters:
r_data (NDArray[np.float64]) – Radial distance from the lidar to the range gate center
signal_strength (NDArray[np.float64]) – Signal strength of the lidar signal (e.g., CNR) at the radial positions r_data
- Returns:
Parameters from the fit: [Rmax, CNRmax, FWHM]
- Return type:
NDArray[np.float64]
- class lidalign.SSC.PeakPulse(gate_length: float = 50)[source]
Bases:
PulseShapePeak pulse shape (Peak to 1 at dr=0). Not realistic, only used for demonstration purposes.
- class lidalign.SSC.PulseShape[source]
Bases:
objectBase Class for pulse shapes of pulses used for wind lidar. Here, the pulse shape means the Range Weighting Function.
Must have a “get_weighting” method implemented in subclass.
- get_inverse_cdf(dr: ndarray[tuple[Any, ...], dtype[float64]]) ndarray[tuple[Any, ...], dtype[float64]][source]
get the inverse cummulative distribution function (CDF) of a pulse shape.
- Parameters:
dr (NDArray[float]) – radial distance to the pulse/range gate center. Resolution of r should be high enough for proper integration. The weighting function must be centered around 0
- Returns:
CDF weighting function (discretized)
- Return type:
NDArray[float]
- class lidalign.SSC.SSC(ds: Dataset, verbose: int = 1)[source]
Bases:
objectBase class for the SeaSurfaceCalibration from surrounding lidar scans and derived lidar-water ranges. For a detailed description of the SSC see the documentation or:
Meyer, P., Rott, A., Schneemann, J., Gramitzky, K., Pauscher, L., & Kühn, M. (2026). Experimental validation of the Sea-Surface-Calibration for scanning lidar static elevation offset determination (in preparation).
Many parts are based on the Sea Surface Leveling by Rott et al. 2022: Rott, A., Schneemann, J., Theuer, F., Trujillo Quintero, J. J., & Kühn, M. (2022). Alignment of scanning lidars in offshore wind farms. Wind Energy Science, 7(1), 283–297. https://doi.org/10.5194/wes-7-283-2022
- get_all_water_ranges(n_processes: int = 1, pulse: PulseShape = None, saveplot: bool | str = False, **kwargs)[source]
Wrapper for WaterRangeDetermination objects to get the lidar-water range for all measurements in the dataset.
- Parameters:
n_processes (int, optional) – number of processes to use for parallel processing. Defaults to 1.
pulse (PulseShape, optional) – pulse shape object. Only used for method LinSig. Defaults to None.
saveplot (bool | str, optional) – whether to save plots or the path to save them. Defaults to False.
**kwargs – additional arguments to pass to WaterRangeDetection.get_water_range_from_cnr()
- Returns:
SSC object with updated distance_ds attribute containing the lidar-water ranges
- Return type:
self
- static get_misalignment(data_in: Dataset, consider_elevation_offset: bool = True, plot: bool = False, print_help: bool = True, fit_method: Literal['lorentz', 'LSQ'] = 'lorentz', consider_earth_curvature: bool = True, error_normalized: bool = True, x0=[0, 0, 20], ax=None, fixed_height: float | None = None, reduce_errors: Literal['ranges', 'elevation'] = 'ranges', return_fit: bool = False) SSCFitResults[source]
Get the external misalignment (pitch, roll, height), and internal misalignment (elevation offset) from multi-elevation scans of the water
Following: Meyer, P., Rott, A., Schneemann, J., Gramitzky, K., Pauscher, L., & Kühn, M. (2026). Experimental validation of the Sea-Surface-Calibration for scanning lidar static elevation offset determination (in preparation).
- Parameters:
data_in (xr.Dataset) – Dataset containing the variables “water_range”, “azimuth”, “elevation” and the dimensions “time”
consider_elevation_offset (bool, optional) – Consider elevation offset. Defaults to True.
plot (bool, optional) – Flag, whether to plot the results. Defaults to False.
print_help (bool, optional) – Flag, whether to print help information about necessary rotation of the legs. Defaults to True.
fit_method (Literal["lorentz", "LSQ"], optional) – Residual summation function to use. Defaults to ‘lorentz’.
consider_earth_curvature (bool, optional) – Consider earth curvature for the sea surface. Defaults to True.
error_normalized (bool, optional) – Normalize errors with range. Defaults to True.
x0 (list, optional) – Initial guess for parameters. Defaults to [0, 0, 20].
ax (_type_, optional) – axes for plotting, if plot=True. Defaults to None.
fixed_height (float | None, optional) – Fixed height value. If this value is set, it will force the fit to use this value. Can be used, if the exact height is known. Otherwise, it is neglected. Defaults to None.
reduce_errors (Literal['ranges','elevation'], optional) – How to calculate the residuals. Defaults to ‘ranges’. Possibly, in Gramitzky 2025 et al. only elevation residuals are used (to be confirmed).
return_fit (bool, optional) – Whether to return the fit data. Defaults to False.
- Returns:
Results of the misalignment fit. Contains the attributes sucess (bool), x (NDArray[np.float64]) with the fitted parameters and result_dict (dict) with the named parameters.
- Return type:
- static interprete_results(result, gewinde_steigung: float = 0.00175, distance_feet: float = 1) None[source]
Interprete the fit results so you can adjust your legs accordingly. Prints the number of rotations for each leg.
- Parameters:
result (dict) – Dictionary containing the fit results with keys “roll”, “pitch”
gewinde_steigung (float, optional) – thread pitch of the lidar legs [m/rotation]. Typically 1.75 mm per rotation. Defaults to 1.75/1000 m/rotation.
distance_feet (float, optional) – distance between the legs of the lidar that contain the screws for lifting [m]. Defaults to 1,.
- static rotated_water_elevation(los_water_range: ndarray[tuple[Any, ...], dtype[float64]], los_azimuth: ndarray[tuple[Any, ...], dtype[float64]], lidar_roll: float, lidar_pitch: float, lidar_height: float, elevation_offset: float = 0) ndarray[tuple[Any, ...], dtype[float64]][source]
Calculate los_elevation with given water range, azimtuh, lidar roll, pitch and height for the simplified and generalized SSL. not used in SSC but in generalized SSL?
Following Gramitzky et al. 2025: Gramitzky, K., Jäger, F., Callies, D., Hildebrand, T., Lundquist, J. K., & Pauscher, L. (2025). Alignment of Scanning Lidars in Offshore Campaigns – an Extension of the Sea Surface Levelling Method. Wind and the atmosphere/Wind and turbulence. https://doi.org/10.5194/wes-2025-191
- Parameters:
los_water_range (NDArray[np.float]) – Lidar-water range in line-of-sight direction [m]
los_azimuth (NDArray[np.float]) – Line-of-sight azimuth angle [deg]
lidar_roll (float) – Lidar roll angle [deg]
lidar_pitch (float) – Lidar pitch angle [deg]
lidar_height (float) – Lidar height above sea surface [m]
elevation_offset (float, optional) – Elevation offset in line-of-sight direction [deg]. Defaults to 0.
- Returns:
Calculated los_elevation of line-of-sight [deg
- Return type:
NDArray[np.float]
- static rotated_water_range(los_elevation: ndarray[tuple[Any, ...], dtype[float64]], los_azimuth: ndarray[tuple[Any, ...], dtype[float64]], lidar_roll: float, lidar_pitch: float, lidar_height: float, los_elevation_offset: float = 0, consider_earth_curvature: bool = True) ndarray[tuple[Any, ...], dtype[float64]][source]
Calculate the lidar-water range for given commanded LOS elevation and azimuth and external misalignment (roll, pitch, lidar_height) and internal misalignment (los_elevation_offset) for the Sea Surface Scans
Following Rott et al. 2022 and Meyer et al. 2026 (their nomenclature is used here)
- Parameters:
los_elevation (NDArray[np.float]) – _line-of-sight elevation angles [deg]
los_azimuth (NDArray[np.float]) – _line-of-sight azimuth angles [deg]
lidar_roll (float) – roll of lidar [deg]. Right handed coordinate system. roll around y-axis (lidar north)
lidar_pitch (float) – pitch of lidar [deg]. Right handed coordinate system. pitch around x-axis (lidar east)
lidar_height (float) – height of lidar above sea surface [m]
los_elevation_offset (float, optional) – internal offset of line-of-sight elevation [deg]. Defaults to 0.
consider_earth_curvature (bool, optional) – Consider earth curvature or flat plate of earth (False). Defaults to True.
- Returns:
lidar-water range [m]
- Return type:
NDArray[np.float]
- class lidalign.SSC.SSCFitResults(success: bool, x: ndarray[tuple[Any, ...], dtype[float64]], result_dict: dict = None, residual: float = inf)[source]
Bases:
objectDataclass for SSC fit results.
- Parameters:
success – bool Indicates if the fitting process was successful.
x – NDArray[np.float64] Array of fitted parameter values.
result_dict – dict, optional Dictionary containing additional fit results and metadata. Default is None.
fig (optional) – Matplotlib figure object containing the fit plot. Default is None.
- residual: float = inf
- result_dict: dict = None
- success: bool
- x: ndarray[tuple[Any, ...], dtype[float64]]
- class lidalign.SSC.WaterRangeDetection(data_db: Dataset, pulse: PulseShape = None, verbose: int = 0, input_in_db: bool = True)[source]
Bases:
objectBase class for the detection of the lidar-water range r_w
- static convolution_fit_error(params: tuple, data: DataArray, pulse: PulseShape, use_linear_scale: bool) float[source]
Calculate the fit error (sum of sqared errors) for convolution with fixed FWHM of pulse
- Parameters:
params (tuple) – Parameters for the fit error calculation (R_water, cnr_slope, cnr_offset, cnr_noise, FWHM)
data (xr.DataArray) – input data for fit error calculation, must have variables “range” and “cnr”
pulse (PulseShape) – pulse shape object
use_linear_scale (bool) – whether the input data is in linear scale
- Returns:
Calculated fit error
- Return type:
float
- static convolution_fit_error_pulsevar(params: tuple, data: DataArray, pulse: PulseShape, use_linear_scale: bool) float[source]
Calculate the fit error betwenen input data and modeled data: model from pulse with variable FWHM
- Parameters:
params (tuple) – Parameters for the fit error calculation (R_water, cnr_slope, cnr_offset, cnr_noise, FWHM)
data (xr.DataArray) – input data for fit error calculation, must have variables “range” and “cnr”
pulse (PulseShape) – pulse shape object
use_linear_scale (bool) – whether the input data is in linear scale
- Returns:
Calculated fit error
- Return type:
float
- get_water_range_from_cnr(min_cnr: float = -22, cnr_hard_target: float = 0, use_bounds: bool = True, func: Literal['LinSig', 'dBSig', 'Gra24', 'Rot22', 'dBConvo', 'Convo_pulse', 'Convo', 'LinConvo_pulse'] = 'LinSig', fit_method: Literal['curve_fit', 'minimize', 'LSQ'] = 'LSQ', dist_guess: float = None, high_cnr_bounds=[-23, -3], low_cnr_bounds=[-40, -20], cnr_noise_cut=-40, growth_rate_bounds=[0.0005, 1], lin_factor_bounds=[-0.015, 0.03], lin_fac_guess: float = 0.01, cnr_noise_first_guess: float = -32, cnr_noise_bounds: list = [-40, -25], _first_guess_roll_window: float = 10, _first_guess_scale: Literal['dB', 'Lin'] = 'Lin', show_plot: bool = False, return_fit: bool = False) CNRFitResult[source]
Get the lidar-water range from the CNR signal of a pulsed lidar.
All Inputs/Boundaries are taken in dB scale.
- Should be called like this:
fit_result = WaterRangeDetection(data_db).get_water_range_from_cnr(…) distance = fit_result.r_water
- Follows:
Meyer, P.J., Rott, A., Schneemann, J., Gramitzky, K., Pauscher, L., & Kühn, M. (2026). Experimental validation of the Sea-Surface-Calibration for scanning lidar static elevation offset determination (in preparation).
- Parameters:
min_cnr (float, optional) – minimum CNR to be considered as valid signal [dB]. If max(CNR) is below, the distance is NAN. Defaults to -22 dB.
cnr_hard_target (float, optional) – CNR value of hard target [dB]. If max(CNR) is larger, hard target is detected and distance is NAN. Defaults to 0 dB.
use_bounds (bool, optional) – Whether to use bounds for the fit parameters. Defaults to True.
func (Literal["LinSig","dBSig", "Gra24", "Rot21","Convo", "Convo_pulse",'LinConvo','LinConvo_pulse'], optional) – Method to use for the fitting. A full description is provided in the documentation. Defaults to “LinSig”.
fit_method (Literal["curve_fit", "minimize",'LSQ'], optional) – Method to fit the least squares. Defaults to “LSQ”.
show_plot (bool, optional) – Show plot of fit afterwards. Defaults to False.
return_fit (bool, optional) – Whether to return the fit result object. Defaults to False.
dist_guess (float, optional) – Initial guess for the distance to support fit [m]. If None, not considered. Defaults to None.
high_cnr_bounds (list, optional) – Bounds for high CNR values in fit [dB]. Defaults to [-23, -3].
low_cnr_bounds (list, optional) – Bounds for low CNR values in fit [dB]. Defaults to [-40, -20].
cnr_noise_cut (int, optional) – Threshold to cut off CNR noise in fit [dB]. Only used for LinSig. Should be adjusted to the actual data. Defaults to -40 dB.
growth_rate_bounds (list, optional) – Bounds for the growth rate parameter of the sigmoid. Defaults to [0.0005,1].
lin_factor_bounds (list, optional) – Bounds for the linear factor parameter of the CNR decay before the water [dB/m]. Defaults to [-0.015, 0.03].
lin_fac_guess (float, optional) – Initial guess for the linear factor [dB/m]. Defaults to 0.01.
cnr_noise_first_guess (float, optional) – Initial guess for the CNR noise level [dB]. Defaults to -32.
cnr_noise_bounds (list, optional) – Bounds for the CNR noise level [dB]. Defaults to [-40,-25].
_first_guess_roll_window (float, optional) – Number of values for rolling mean in CNR to find first guess for middle range, where change is the strongest. This parameter should be chosen accordinly to the spatial resoultion of the range gate centers and the probe volume length. In Principle, the value should be L_probevolume/ dr.
_first_guess_scale – Literal[‘dB’,’Lin’] = ‘Lin’: Whether to use the CNR in dB scale or linear scale for the first guess estimation. Especially for far ranges (low elevation angles) LinSig can lead to a faulty first guess, as the CNR (in linear scale) might drop faster at shorter ranges (due to the atmospheric backscatter). The DB scale however can have strugles, when there is a lot of scatter in the CNR noise -> Should only be used with CNR noise removal parameter set. Defaults to ‘Lin’.
- Returns:
Fit Results dataclass with water range and fit success (and plot or fit data if requested)
- Return type:
- static model_cnr_signal_convolution(r: ndarray, R_water: float, cnr_slope: float, cnr_offset: float, pulse_object: PulseShape, cnr_noise: float = -30, return_dB: bool = False)[source]
Model the cnr signal convolution of a linear decay with a drop at R_water and a pulse shape, given by a pulse_object –> the returned signal is given in the linear scale!
- Parameters:
r (np.ndarray) – Range array where convolution is calculated at - must be equidistant [m]
R_water (float) – Range to water surface [m]
cnr_slope (float) – Slope of the linear decay [dB/m]
cnr_offset (float) – Offset of the linear decay [dB]
pulse_object (PulseShape) – Pulse shape object
cnr_noise (float, optional) – Noise level in CNR. Defaults to -28. [dB]
- Returns:
Modeled CNR signal at range r
- Return type:
NDArray[np.float64]
- lidalign.SSC.db2linear(dbsignal: float | ndarray[tuple[Any, ...], dtype[float64]]) float | ndarray[tuple[Any, ...], dtype[float64]][source]
Convert signal in dB scale into linear scale: S_{lin} = 10^{S_{dB}/10},
- Parameters:
dbsignal (float | NDArray[np.float64]) – signal strength in dB scale
- Returns:
signal strength in linear scale
- Return type:
float|NDArray[np.float64]
- lidalign.SSC.distance_to_water(data, min_cnr, show_plot=0, high_cnr_ub=0, high_cnr_lb=-25, low_cnr_ub=-15, low_cnr_lb=-30, distance_ub=3000, distance_lb=200, growth_ub=1, growth_lb=0)[source]
- lidalign.SSC.inverse_sigmoid(r: float | ndarray[tuple[Any, ...], dtype[float64]], mid: float, down: float, up: float, growth: float) float | ndarray[tuple[Any, ...], dtype[float64]][source]
Sigmoid function for fitting CNR data for water range detection.
Following: Rott, A., Schneemann, J., Theuer, F., Trujillo Quintero, J. J., & Kühn, M. (2022). Alignment of scanning lidars in offshore wind farms. Wind Energy Science, 7(1), 283–297. https://doi.org/10.5194/wes-7-283-2022
- Parameters:
r (float|NDArray[np.float64]) – radial distance of range gate center
up (float) – upper CNR level before drop [dB]
down (float) – lower CNR level after drop [dB]
mid (float) – radial distance of the inflection point of the sigmoid -> corresponds to water range [m]
growth (float) – growth rate of the sigmoid function [-]
- Returns:
Sigmoid function (r)
- lidalign.SSC.inverse_sigmoid_Gra24(r: float | ndarray[tuple[Any, ...], dtype[float64]], mid: float, down: float, up: float, growth: float, lin_fac: float) ndarray[tuple[Any, ...], dtype[float64]][source]
Sigmoid function for fitting CNR data for water range detection, using a linear decay/increase before the drop.
In contrast to the original function from Rott et al. (2022), a linear factor is added to account for a linear decrease before the drop.
Following: Gramitzky, K., Jäger, F., Hildebrand, T., Gloria, N., Riechert, J., Steger, M., & Pauscher, L. (2024). Alignment calibration and correction for offshore wind measurements using scanning lidars. Journal of Physics: Conference Series, 2767(4), 042014. https://doi.org/10.1088/1742-6596/2767/4/042014
- Parameters:
r (float|NDArray[np.float64]) – radial distance of range gate center
mid (float) – radial distance of the inflection point of the sigmoid -> corresponds to water range [m]
down (float) – lower CNR level after drop [dB]
up (float) – upper CNR level before drop [dB]
growth (float) – growth rate of the sigmoid function [-]
lin_fac (float) – linear factor to adjust linear decrease before the drop [dB/m]
- Returns:
sigmoid value, depending on input type of position
- lidalign.SSC.inverse_sigmoid_dbscale(r: float | ndarray[tuple[Any, ...], dtype[float64]], mid: float, down: float, up: float, growth: float, lin_fac: float) ndarray[tuple[Any, ...], dtype[float64]][source]
Inverse sigmoid in DB scale. All inputs must be in dB scale.
we do not recommend using this function
Following: Meyer, P., Rott, A., Schneemann, J., Gramitzky, K., Pauscher, L., & Kühn, M. (2026). Experimental validation of the Sea-Surface-Calibration for scanning lidar static elevation offset determination (in preparation).
- Parameters:
r (float | NDArray[np.float64]) – radial distance of range gate center
up (float) – upper CNR level before drop [dB]
down (float) – lower CNR level after drop [dB]
mid (float) – radial distance of the inflection point of the sigmoid -> corresponds to water range [m]
growth (float) – growth rate of the sigmoid function [-]
lin_fac (float) – linear factor to adjust linear decrease before the drop [dB/m]
- Returns:
Sigmoid function (r) in dB scale
- Return type:
NDArray[np.float64]
- lidalign.SSC.inverse_sigmoid_linscale(r: float | ndarray[tuple[Any, ...], dtype[float64]], mid: float, down: float, up: float, growth: float, lin_fac: float) ndarray[tuple[Any, ...], dtype[float64]][source]
Inverse sigmoid in linear scale. Same as inverse_sigmoid_dB but output format is linear scale. INPUT MUST BE dB SCALE
- Following:
Meyer, P., Rott, A., Schneemann, J., Gramitzky, K., Pauscher, L., & Kühn, M. (2026). Experimental validation of the Sea-Surface-Calibration for scanning lidar static elevation offset determination (in preparation).
- Parameters:
r (float | NDArray[np.float64]) – radial distance of range gate center
up (float) – upper CNR level before drop [dB]
down (float) – lower CNR level after drop [dB]
mid (float) – radial distance of the inflection point of the sigmoid -> corresponds to water range [m]
growth (float) – growth rate of the sigmoid function [-]
lin_fac (float) – linear factor to adjust linear decrease before the drop [dB/m]
- Returns:
Sigmoid function (r) in linear scale
- Return type:
NDArray[np.float64]
- lidalign.SSC.linear2db(linearsignal: float | ndarray[tuple[Any, ...], dtype[float64]]) float | ndarray[tuple[Any, ...], dtype[float64]][source]
Convert signal in linear scale into dB scale: S_{dB} = 10 log_{10}(S_{lin}),
- Parameters:
linearsignal (float | NDArray[np.float64]) – signal strength in linear scale
- Returns:
signal strength in dB scale
- Return type:
float|NDArray[np.float64]
- lidalign.SSC.model_cnr_signal_CDF(r, R_water, uppervalue: float, lowervalue: float, lin_fac: float, pulse_object: PulseShape)[source]
Model the CNR signal using the inverse CDF of a pulse shape. Not really used anymore
lidalign.hard_target_elevation_mapping module
- class lidalign.hard_target_elevation_mapping.HardTargetMappingElevation(azimuth: ndarray, delta_elevation: ndarray, unc_azimuth: ndarray | None = None, unc_elevation: ndarray | None = None)[source]
Bases:
object- fit(n_mc=1000, typ: Literal['cosine', 'pitchroll'] = 'pitchroll')[source]
Fit the data to the cosine curve.
- class lidalign.hard_target_elevation_mapping.MonteCarloFunc(parameters: list[DataFrame], n: int = 1000)[source]
Bases:
object
- lidalign.hard_target_elevation_mapping.cosine_curve(azi: float | ndarray, phase_offset: float, amplitude: float, offset: float) float | ndarray[source]
Cosine curve
- Parameters:
azi (float | np.ndarray) – Azimuth angle [degree]
phase_offset (float) – Phase Offset of cosine curve [degree]
amplitude (float) – Amplitude of cosine curve
offset (float) – Y-offset of cosine curve
- Returns:
Y values of cosine curve
- Return type:
float|np.ndarray
- lidalign.hard_target_elevation_mapping.uncertain_df(mean_values: ndarray, uncertainty_params: ndarray, uncertainty_distribution: str = 'normal') DataFrame[source]
Define a Uncertainty dataframe, that contains the required statistics
- Parameters:
mean_values (np.ndarray) – mean values of a variable
uncertainty_params (np.ndarray) – uncertainty parameters, here often standard uncertainty of a variable
uncertainty_distribution (str, optional) – How is the uncertainty distributed?. Defaults to “normal”/gaussian.
- Returns:
Dataframe containing the Uncertainty informations
- Return type:
pd.DataFrame
lidalign.io module
Scripts related to reading of WindCubeScan data from Vaisala. Authors: Paul Meyer, Ignace Ransquin (ForWind WESys) Created: 2025-08-25
V0 2025-08-25: Initial push to Git + minor changes
- class lidalign.io.FileDB(path: str, regex: str = '*.*')[source]
Bases:
object- filter_file_names(start: str = None, end: str = None, timedelta_back: str = None, closest_to_time: str = None, max_timedist='3h')[source]
- get_filtered_filelist(start: str = None, end: str = None, timedelta_back: str = None, filename_regex: str = None)[source]
Get the reduced filelist, depending on start and end time, as well as a regex to search for in the filenames
- Parameters:
start (str, optional) – Start date+time of the lidar data to read. If None, start is set to first date. Defaults to None.
end (str, optional) – End date+time of lidar data to read. If None, end is set to NOW. Defaults to None.
filename_regex (str, optional) – regex to filter the filenames. Can be e.g., the scan type and number: ppi_93 or similar. Defaults to None.
- Returns:
Database Object with filtered_files_list attribute, which contains the filtered files list
- Return type:
self
- class lidalign.io.RawEnvironmentalDB(path: str)[source]
Bases:
FileDB- filepattern = '.*(\\d{4}\\d{2}\\d{2}_\\d{2}\\d{2}\\d{2}).csv'
- names = {'Time': <class 'str'>}
- class lidalign.io.WindCubeScanDB(path_to_campaign: str, datatype: Literal['wind_and_aerosols_data', 'environmental_data', 'scans'] = 'wind_and_aerosols_data', file_structure: Literal['native_vaisala', 'flat'] = 'native_vaisala', verbose=0, prefilter_dates: dict | None = None, position: list[float, float, float] = None)[source]
Bases:
FileDBSetup of “database” to read measurements from Vaisala WindCube Scan devices. Author: Paul Meyer, Ignace Ransquin (WESys) Created: 2025-08-25
- get_data(start: str = None, end: str = None, concatenated=False, timedelta_back: str = None, filename_regex: str = None, **read_kwargs)[source]
Wrapper to get the data after filtering the total files list
- Parameters:
concatenated (bool, optional) – If true, retrieved datasets are concatenated to a single one with time as main dimension. If False, a list of datasets is returned. Defaults to False.
- kwargs:
keywords for WindCubeScanDB.get_filtered_filelist()
- Returns:
List or Xarray Dataset, depending on concatenated keyword
- Return type:
xr.dataset or List
- read_wind_files(files: list, concatenated=True, generate_utc_time=True, sel_kwargs: {} = None, query_kwargs: {} = None, max_n: int = None, processes: int = 1, VAD_ranges: list = None, **read_kwargs)[source]
Read datasets from filenames
- Parameters:
files (list) – list of files to read
concatenated (bool, optional) – If true, retrieved datasets are concatenated to a single one with time as main dimension. If False, a list of datasets is returned. Defaults to False.
generate_utc_time (bool, optional) – If true, an additional variable with pandas datetime is created. Defaults to True.
- Returns:
depending on concatenated, list or concatenated dataset of lidar measuremnts
- Return type:
data
lidalign.north_alignment module
The scripts here are based on the scripts by Rott et al 2022: Rott, A., Schneemann, J., Theuer, F., Trujillo Quintero, J. J., & Kühn, M. (2022). Alignment of scanning lidars in offshore wind farms. Wind Energy Science, 7(1), 283–297. https://doi.org/10.5194/wes-7-283-2022
They are also available at: Andreas Rott, Jörge Schneemann, & Frauke Theuer. (2021). AndreasRott/Alignment_of_scanning_lidars_in_offshore_wind_farms: Version1.0 (Release1.0.0). Zenodo. https://doi.org/10.5281/zenodo.5654919
- class lidalign.north_alignment.Northalignment(HardTargets_coordinates: DataFrame)[source]
Bases:
objectClass for the north alignment of a scanning lidar based on surrounding hard target measurements.
- static cost_function(params: tuple[float, float, float], HardTargets: DataFrame, lidardata: DataFrame, summed_cost: bool = True)[source]
Cost function to optimize the lidar positioning and alignment.
- Reference:
Rott, A., Schneemann, J., Theuer, F., Trujillo Quintero, J. J., & Kühn, M. (2022). Alignment of scanning lidars in offshore wind farms. Wind Energy Science, 7(1), 283–297. https://doi.org/10.5194/wes-7-283-2022
- Parameters:
params (list[float, float, float]) – Includes x_offset, y_offset and azi_offset (each float) of the lidar in the global coordinate system.
HardTargets (pd.DataFrame) – DataFrame containing the positions of the Hard Targets. Columns must be x and y
lidardata (pd.DataFrame) – lidar data from the hard target measurements. Must contain columns range and azimuth. Here, no elevation is considered.
- Returns:
cost of the function with given parameters
- Return type:
float
- fit(lidardata: Dataset, initial_guess: tuple[float, float, float], CNR_hardtarget: float = 0, method_use: Literal['V1', 'V2'] = 'V2', R_add_hardtargets: float = 0, bounds_dist: float = 300, bounds_angle: float = 30, redo_iterations: int = 2, max_distance: float = 100, print_result: bool = False, plot: bool = False)[source]
Fit the measurements with the coordinate list
- Parameters:
lidardata (xr.Dataset) – xarray dataset with the columns cnr [db], azimuth [°] and range [m], containing the measurements of the hard targets. Here, no elevation is considered.
initial_guess (tuple[float, float, float]) – Initial guess [position x, position y, initial lidar orientation] for the fit.
CNR_hardtarget (float, optional) – CNR threshold for a measurement to be considered a hard target measurement, in dB. Defaults to 0.
method_use (Literal["V1", "V2"], optional) – Method to use for the fit, V1 is Rott et al. 2022, V2 is the updated version. Defaults to “V2”.
R_add_hardtargets (float, optional) – Range to add to the detected hard target positions, as the turbine coordinate is for the center of the tower, while the hard target results in the tower shell position. Defaults to 5.
bounds_dist (float, optional) – bounds for the posotion to deviate from the initial guess. Defaults to 300.
bounds_angle (float, optional) – bounds for the posotion to deviate from the initial guess. Defaults to 30.
redo_iterations (int, optional) – number of iterations to redo the fit with outliers removed. Defaults to 2.
max_distance (float, optional) – maximum distance for a measurement to be considered an outlier. Defaults to 100.
print_result (bool, optional) – whether to print the fit result. Defaults to False.
plot (bool, optional) – whether to plot the fit result. Defaults to False.
Returns – result: result of the fit, including the optimized parameters and the cost function value. Object of scipy fit results
- static plot(params: tuple[float, float, float], lidar_data_df: DataFrame, HardTargets: DataFrame, interactive: bool = False, all_hard_targets_lidar: DataFrame | None = None)[source]
Plot the alignment and results
- Parameters:
params (tuple[float, float, float]) – params of the fit, including x_offset, y_offset and azi_offset (each float) of the lidar in the global coordinate system.
lidar_data_df (pd.DataFrame) – lidar data, must contain columns range and azimuth. Here, no elevation is considered.
HardTargets (pd.DataFrame) – Hard target positions, must contain columns x and y
interactive (bool, optional) – whether to create an interactive plot (with plotly). Defaults to False.
all_hard_targets_lidar (pd.DataFrame | None, optional) – all detected hard targets. Defaults to None.
lidalign.utils module
- lidalign.utils.publication_figure(relative_width: float = 0.5, width_paper: float = 6.30045, height=3, fig_only: bool = False, **kwargs)[source]
width_paper = 6.3.. -> width of IOP template, find more here: https://tex.stackexchange.com/questions/39383/determine-text-width