Loads Module
The loads package of the MHKiT (Marine and Hydrokinetic Toolkit) library provides tools and functionalities for analyzing and visualizing loads data from marine and hydrokinetic (MHK) devices. This package is designed to assist engineers, researchers, and analysts in understanding the forces and stresses applied to MHK devices under various operational and environmental conditions.
General
This module provides tools for analyzing and processing data signals related to turbine blade performance and fatigue analysis. It implements methodologies based on standards such as IEC TS 62600-3:2020 ED1, incorporating statistical binning, moment calculations, and fatigue damage estimation using the rainflow counting algorithm. Key functionalities include:
bin_statistics: Bins time-series data against a specified signal, such as wind speed, to calculate mean and standard deviation statistics for each bin, following IEC TS 62600-3:2020 ED1 guidelines. It supports output in both pandas DataFrame and xarray Dataset formats.
blade_moments: Calculates the flapwise and edgewise moments of turbine blades using derived calibration coefficients and raw strain signals. This function is crucial for understanding the loading and performance characteristics of turbine blades.
damage_equivalent_load: Estimates the damage equivalent load (DEL) of a single data signal using a 4-point rainflow counting algorithm. This method is vital for assessing fatigue life and durability of materials under variable amplitude loading.
References: - C. Amzallag et. al., International Journal of Fatigue, 16 (1994) 287-293. - ISO 12110-2, Metallic materials - Fatigue testing - Variable amplitude fatigue testing. - G. Marsh et. al., International Journal of Fatigue, 82 (2016) 757-765.
Bins calculated statistics against data signal (or channel) according to IEC TS 62600-3:2020 ED1. |
|
Transfer function for deriving blade flap and edge moments using blade matrix. |
|
Calculates the damage equivalent load of a single data signal (or channel) based on IEC TS 62600-3:2020 ED1. |
- mhkit.loads.general.bin_statistics(data: DataFrame | Dataset, bin_against: ndarray, bin_edges: ndarray, data_signal: List[str] | None = None, to_pandas: bool = True) Tuple[DataFrame | Dataset, DataFrame | Dataset] [source]
Bins calculated statistics against data signal (or channel) according to IEC TS 62600-3:2020 ED1.
- Parameters:
data (pandas DataFrame or xarray Dataset) – Time-series statistics of data signal(s)
bin_against (array) – Data signal to bin data against (e.g. wind speed)
bin_edges (array) – Bin edges with consistent step size
data_signal (list, optional) – List of data signal(s) to bin, default = all data signals
to_pandas (bool (optional)) – Flag to output pandas instead of xarray. Default = True.
- Returns:
bin_mean (pandas DataFrame or xarray Dataset) – Mean of each bin
bin_std (pandas DataFrame or xarray Dataset) – Standard deviation of each bim
- mhkit.loads.general.blade_moments(blade_coefficients: ndarray, flap_offset: float, flap_raw: ndarray, edge_offset: float, edge_raw: ndarray) Tuple[ndarray, ndarray] [source]
Transfer function for deriving blade flap and edge moments using blade matrix.
- Parameters:
blade_coefficients (numpy array) – Derived blade calibration coefficients listed in order of D1, D2, D3, D4
flap_offset (float) – Derived offset of raw flap signal obtained during calibration process
flap_raw (numpy array) – Raw strain signal of blade in the flapwise direction
edge_offset (float) – Derived offset of raw edge signal obtained during calibration process
edge_raw (numpy array) – Raw strain signal of blade in the edgewise direction
- Returns:
M_flap (numpy array) – Blade flapwise moment in SI units
M_edge (numpy array) – Blade edgewise moment in SI units
- mhkit.loads.general.damage_equivalent_load(data_signal: ndarray, m: float | int, bin_num: int = 100, data_length: float | int = 600) float [source]
Calculates the damage equivalent load of a single data signal (or channel) based on IEC TS 62600-3:2020 ED1. 4-point rainflow counting algorithm from fatpack module is based on the following resources:
C. Amzallag et. al. Standardization of the rainflow counting method for fatigue analysis. International Journal of Fatigue, 16 (1994) 287-293
ISO 12110-2, Metallic materials - Fatigue testing - Variable amplitude fatigue testing.
G. Marsh et. al. Review and application of Rainflow residue processing techniques for accurate fatigue damage estimation. International Journal of Fatigue, 82 (2016) 757-765
Parameters:
- data_signalarray
Data signal being analyzed
- mfloat/int
Fatigue slope factor of material
- bin_numint
Number of bins for rainflow counting method (minimum=100)
- data_lengthfloat/int
Length of measured data (seconds)
- returns:
DEL (float) – Damage equivalent load (DEL) of single data signal
Graphics
This module provides functionalities for plotting statistical data related to a given variable or dataset.
plot_statistics is designed to plot raw statistical measures (mean, maximum, minimum, and optional standard deviation) of a variable across a series of x-axis values. It allows for customization of plot labels, title, and saving the plot to a file.
plot_bin_statistics extends these capabilities to binned data, offering a way to visualize binned statistics (mean, maximum, minimum) along with their respective standard deviations. This function also supports label and title customization, as well as saving the plot to a specified path.
Plot showing standard raw statistics of variable |
|
Plot showing standard binned statistics of single variable |
- mhkit.loads.graphics.plot_statistics(x: ndarray, y_mean: ndarray, y_max: ndarray, y_min: ndarray, y_stdev: ndarray | None = None, **kwargs: Dict[str, Any]) Axes [source]
Plot showing standard raw statistics of variable
- Parameters:
x (numpy array) – Array of x-axis values
y_mean (numpy array) – Array of mean statistical values of variable
y_max (numpy array) – Array of max statistical values of variable
y_min (numpy array) – Array of min statistical values of variable
y_stdev (numpy array, optional) – Array of standard deviation statistical values of variable
**kwargs (optional) –
- x_labelstring
x axis label for plot
- y_labelstring
y axis label for plot
- titlestring, optional
Title for plot
- save_pathstring
Path and filename to save figure.
- Returns:
ax (matplotlib pyplot axes)
- mhkit.loads.graphics.plot_bin_statistics(bin_centers: ndarray, bin_mean: ndarray, bin_max: ndarray, bin_min: ndarray, bin_mean_std: ndarray, bin_max_std: ndarray, bin_min_std: ndarray, **kwargs: Dict[str, Any]) Axes [source]
Plot showing standard binned statistics of single variable
- Parameters:
bin_centers (numpy array) – x-axis bin center values
bin_mean (numpy array) – Binned mean statistical values of variable
bin_max (numpy array) – Binned max statistical values of variable
bin_min (numpy array) – Binned min statistical values of variable
bin_mean_std (numpy array) – Standard deviations of mean binned statistics
bin_max_std (numpy array) – Standard deviations of max binned statistics
bin_min_std (numpy array) – Standard deviations of min binned statistics
**kwargs (optional) –
- x_labelstring
x axis label for plot
- y_labelstring
y axis label for plot
- titlestring, optional
Title for plot
- save_pathstring
Path and filename to save figure.
- Returns:
ax (matplotlib pyplot axes)
Extreme
This package provides tools and functions for extreme value analysis and wave data statistics.
It includes methods for calculating peaks over threshold, estimating short-term extreme distributions,and performing wave amplitude normalization for most likely extreme response analysis.
Estimate the short-term extreme distribution from the peaks distribution. |
|
Find the block maxima of a time-series. |
|
Approximate the short-term extreme distribution using the block maxima method and the Generalized Extreme Value distribution. |
|
Approximate the short-term extreme distribution using the block maxima method and the Gumbel (right) distribution. |
|
Alias for short_term_extreme. |
|
Approximate the short-term extreme distribution from a timeseries of the response using chosen method. |
|
Return the long-term extreme distribution of a response of interest using the full sea state approach. |
|
Calculate MLER (most likely extreme response) coefficients from a sea state spectrum and a response RAO. |
|
Define the simulation parameters that are used in various MLER functionalities. |
|
Function that renormalizes the incoming amplitude of the MLER wave to the desired peak height (peak to MSL). |
|
Generate the wave amplitude time series at X0 from the calculated MLER coefficients |
|
Find the global peaks of a zero-centered response time-series. |
|
Estimate the number of peaks in a specified period. |
|
Estimate the peaks distribution by fitting a Weibull distribution to the peaks of the response. |
|
Estimate the peaks distribution using the Weibull tail fit method. |
|
Find the best significant wave height threshold for the peaks-over-threshold method. |
|
Estimate the peaks distribution using the peaks over threshold method. |
|
Calculate the value from a given distribution corresponding to a particular return year. |
- mhkit.loads.extreme.ste_peaks(peaks_distribution: rv_continuous, npeaks: float) rv_continuous [source]
Estimate the short-term extreme distribution from the peaks distribution.
- Parameters:
peaks_distribution (scipy.stats.rv_frozen) – Probability distribution of the peaks.
npeaks (float) – Number of peaks in short term period.
- Returns:
short_term_extreme (scipy.stats.rv_frozen) – Short-term extreme distribution.
- mhkit.loads.extreme.block_maxima(time: ndarray, global_peaks_data: ndarray, time_st: float) ndarray [source]
Find the block maxima of a time-series.
The timeseries (time, global_peaks) is divided into blocks of length t_st, and the maxima of each bloock is returned.
- Parameters:
time (np.array) – Time array.
global_peaks_data (np.array) – global peaks timeseries.
time_st (float) – Short-term period.
- Returns:
block_max (np.array) – Block maxima (i.e. largest peak in each block).
- mhkit.loads.extreme.ste_block_maxima_gev(block_max)[source]
Approximate the short-term extreme distribution using the block maxima method and the Generalized Extreme Value distribution.
- Parameters:
block_max (np.array) – Block maxima (i.e. largest peak in each block).
- Returns:
short_term_extreme_rv (scipy.stats.rv_frozen) – Short-term extreme distribution.
- mhkit.loads.extreme.ste_block_maxima_gumbel(block_max)[source]
Approximate the short-term extreme distribution using the block maxima method and the Gumbel (right) distribution.
- Parameters:
block_max (np.array) – Block maxima (i.e. largest peak in each block).
- Returns:
ste (scipy.stats.rv_frozen) – Short-term extreme distribution.
- mhkit.loads.extreme.ste(time: ndarray, data: ndarray, t_st: float, method: str) rv_continuous [source]
Alias for short_term_extreme.
- mhkit.loads.extreme.short_term_extreme(time: ndarray, data: ndarray, t_st: float, method: str) rv_continuous | None [source]
Approximate the short-term extreme distribution from a timeseries of the response using chosen method.
The availabe methods are: ‘peaks_weibull’, ‘peaks_weibull_tail_fit’, ‘peaks_over_threshold’, ‘block_maxima_gev’, and ‘block_maxima_gumbel’. For the block maxima methods the timeseries needs to be many times longer than the short-term period. For the peak-fitting methods the timeseries can be of arbitrary length.
- Parameters:
time (np.array) – Time array.
data (np.array) – Response timeseries.
t_st (float) – Short-term period.
method (string) – Method for estimating the short-term extreme distribution.
- Returns:
short_term_extreme_dist (scipy.stats.rv_frozen) – Short-term extreme distribution.
- mhkit.loads.extreme.full_seastate_long_term_extreme(short_term_extreme_dist, weights)[source]
Return the long-term extreme distribution of a response of interest using the full sea state approach.
- Parameters:
ste (list[scipy.stats.rv_frozen]) – Short-term extreme distribution of the quantity of interest for each sample sea state.
weights (list, np.ndarray) – The weights from the full sea state sampling
- Returns:
ste (scipy.stats.rv_frozen) – Short-term extreme distribution.
- mhkit.loads.extreme.mler_coefficients(rao: ndarray[Any, dtype[float64]] | Series | List[float] | List[int] | DataArray, wave_spectrum: Series | DataFrame | DataArray | Dataset, response_desired: int | float, frequency_dimension: str = '', to_pandas: bool = True) DataFrame | Dataset [source]
Calculate MLER (most likely extreme response) coefficients from a sea state spectrum and a response RAO.
- Parameters:
rao (numpy ndarray) – Response amplitude operator.
wave_spectrum (pandas Series, pandas DataFrame, xarray DataArray, or xarray Dataset) – Wave spectral density [m^2/Hz] indexed by frequency [Hz]. DataFrame and Dataset inputs should only have one data variable
response_desired (int or float) – Desired response, units should correspond to a motion RAO or units of force for a force RAO.
frequency_dimension (string (optional)) – Name of the xarray dimension corresponding to frequency. If not supplied, defaults to the first dimension. Does not affect pandas input.
to_pandas (bool (optional)) – Flag to output pandas instead of xarray. Default = True.
- Returns:
mler (pandas DataFrame or xarray Dataset) – DataFrame containing conditioned wave spectral amplitude coefficient [m^2-s], and Phase [rad] indexed by freq [Hz].
- mhkit.loads.extreme.mler_simulation(parameters: Dict[str, float | int | ndarray] | None = None) Dict[str, float | int | ndarray] [source]
Define the simulation parameters that are used in various MLER functionalities.
See extreme_response_contour_example.ipynb example for how this is useful. If no input is given, then default values are returned.
- Parameters:
parameters (dict (optional)) – Simulation parameters. Keys: —– - ‘startTime’: starting time [s] - ‘endTime’: ending time [s] - ‘dT’: time-step size [s] - ‘T0’: time of maximum event [s] - ‘startx’: start of simulation space [m] - ‘endX’: end of simulation space [m] - ‘dX’: horizontal spacing [m] - ‘X’: position of maximum event [m] The following keys are calculated from the above parameters: - ‘maxIT’: int, maximum timestep index - ‘T’: np.ndarray, time array - ‘maxIX’: int, maximum index for space - ‘X’: np.ndarray, space array
- Returns:
sim (dict) – Simulation parameters including spatial and time calculated arrays.
- mhkit.loads.extreme.mler_wave_amp_normalize(wave_amp: float, mler: DataFrame | Dataset, sim: Dict[str, float | int | ndarray], k: ndarray[Any, dtype[float64]] | List[float] | Series, **kwargs: Any) DataFrame | Dataset [source]
Function that renormalizes the incoming amplitude of the MLER wave to the desired peak height (peak to MSL).
- Parameters:
wave_amp (float) – Desired wave amplitude (peak to MSL).
mler (pandas DataFrame or xarray Dataset) – MLER coefficients generated by ‘mler_coefficients’ function.
sim (dict) – Simulation parameters formatted by output from ‘mler_simulation’.
k (numpy ndarray) – Wave number
frequency_dimension (string (optional)) – Name of the xarray dimension corresponding to frequency. If not supplied, defaults to the first dimension. Does not affect pandas input.
to_pandas (bool (optional)) – Flag to output pandas instead of xarray. Default = True.
- Returns:
mler_norm (pandas DataFrame or xarray Dataset) – MLER coefficients
- mhkit.loads.extreme.mler_export_time_series(rao: ndarray[Any, dtype[float64]] | List[float] | Series, mler: DataFrame | Dataset, sim: Dict[str, float | int | ndarray], k: ndarray[Any, dtype[float64]] | List[float] | Series, **kwargs: Any) DataFrame | Dataset [source]
Generate the wave amplitude time series at X0 from the calculated MLER coefficients
- Parameters:
rao (numpy ndarray) – Response amplitude operator.
mler (pandas DataFrame or xarray Dataset) – MLER coefficients dataframe generated from an MLER function.
sim (dict) – Simulation parameters formatted by output from ‘mler_simulation’.
k (numpy ndarray) – Wave number.
frequency_dimension (string (optional)) – Name of the xarray dimension corresponding to frequency. If not supplied, defaults to the first dimension. Does not affect pandas input.
to_pandas (bool (optional)) – Flag to output pandas instead of xarray. Default = True.
- Returns:
mler_ts (pandas DataFrame or xarray Dataset) – Time series of wave height [m] and linear response [*] indexed by time [s].
- mhkit.loads.extreme.global_peaks(time: ndarray, data: ndarray) Tuple[ndarray, ndarray] [source]
Find the global peaks of a zero-centered response time-series.
The global peaks are the maxima between consecutive zero up-crossings.
- Parameters:
time (np.array) – Time array.
data (np.array) – Response time-series.
- Returns:
time_peaks (np.array) – Time array for peaks
peaks (np.array) – Peak values of the response time-series
- mhkit.loads.extreme.number_of_short_term_peaks(n_peaks: int, time: float, time_st: float) float [source]
Estimate the number of peaks in a specified period.
- Parameters:
n_peaks (int) – Number of peaks in analyzed timeseries.
time (float) – Length of time of analyzed timeseries.
time_st (float) – Short-term period for which to estimate the number of peaks.
- Returns:
n_st (float) – Number of peaks in short term period.
- mhkit.loads.extreme.peaks_distribution_weibull(peaks_data: ndarray[Any, dtype[float64]]) rv_continuous [source]
Estimate the peaks distribution by fitting a Weibull distribution to the peaks of the response.
The fitted parameters can be accessed through the params field of the returned distribution.
- Parameters:
peaks_data (NDArray[np.float64]) – Global peaks.
- Returns:
peaks (scipy.stats.rv_frozen) – Probability distribution of the peaks.
- mhkit.loads.extreme.peaks_distribution_weibull_tail_fit(peaks_data: ndarray[Any, dtype[float64]]) rv_continuous [source]
Estimate the peaks distribution using the Weibull tail fit method.
The fitted parameters can be accessed through the params field of the returned distribution.
- Parameters:
peaks_data (np.array) – Global peaks.
- Returns:
peaks (scipy.stats.rv_frozen) – Probability distribution of the peaks.
- mhkit.loads.extreme.automatic_hs_threshold(peaks: ndarray[Any, dtype[float64]], sampling_rate: float, initial_threshold_range: Tuple[float, float, float] = (0.99, 0.995, 0.001), max_refinement: int = 5) Tuple[float, float] [source]
Find the best significant wave height threshold for the peaks-over-threshold method.
This method was developed by:
> Neary, V. S., S. Ahn, B. E. Seng, M. N. Allahdadi, T. Wang, Z. Yang and R. He (2020). > “Characterization of Extreme Wave Conditions for Wave Energy Converter Design and > Project Risk Assessment.” > J. Mar. Sci. Eng. 2020, 8(4), 289; https://doi.org/10.3390/jmse8040289.
Please cite this paper if using this method.
After all thresholds in the initial range are evaluated, the search range is refined around the optimal point until either (i) there is minimal change from the previous refinement results, (ii) the number of data points become smaller than about 1 per year, or (iii) the maximum number of iterations is reached.
- Parameters:
peaks (NDArray[np.float64]) – Peak values of the response time-series.
sampling_rate (float) – Sampling rate in hours.
initial_threshold_range (Tuple[float, float, float]) – Initial range of thresholds to search. Described as (min, max, step).
max_refinement (int) – Maximum number of times to refine the search range.
- Returns:
Tuple[float, float] – The best threshold and its corresponding unit.
- mhkit.loads.extreme.peaks_distribution_peaks_over_threshold(peaks_data: ndarray[Any, dtype[float64]], threshold: float | None = None) rv_continuous [source]
Estimate the peaks distribution using the peaks over threshold method.
This fits a generalized Pareto distribution to all the peaks above the specified threshold. The distribution is only defined for values above the threshold and therefore cannot be used to obtain integral metrics such as the expected value. A typical choice of threshold is 1.4 standard deviations above the mean. The peaks over threshold distribution can be accessed through the pot field of the returned peaks distribution.
- Parameters:
peaks_data (NDArray[np.float64]) – Global peaks.
threshold (Optional[float]) – Threshold value. Only peaks above this value will be used. Default value calculated as: np.mean(x) + 1.4 * np.std(x)
- Returns:
peaks (rv_continuous) – Probability distribution of the peaks.
- mhkit.loads.extreme.return_year_value(ppf: Callable[[float], float], return_year: float, short_term_period_hr: float) float [source]
Calculate the value from a given distribution corresponding to a particular return year.
- Parameters:
ppf (callable function of 1 argument) – Percentage Point Function (inverse CDF) of short term distribution.
return_year (int, float) – Return period in years.
short_term_period_hr (int, float) – Short term period the distribution is created from in hours.
- Returns:
value (float) – The value corresponding to the return period from the distribution.