Wave Module
The wave module contains a set of functions to calculate quantities of interest for wave energy converters (WEC).
Note
The names of the functions below are of the convention path.path.function. Only the function name is used when calling the function in MATLAB. For example, to call on mhkit.wave.io.read_NDBC_file simply use read_NDBC_file.
IO
The io submodule contains the following functions to download and load National Data Buoy Center (NDBC), including real time and historical data, and SWAN model data into structures.
Functions |
Description |
|---|---|
|
Reads a NDBC wave buoy data file (from https://www.ndbc.noaa.gov) into a structure. |
|
Returns the NDBC stations IDs, years, and file names for a requested parameter. |
|
Returns requested NDBC data from passed filenames and parameter. |
|
Reads in SWAN ASCII block format output and returns a data structure. |
|
Reads in SWAN ASCII table format output and returns a data structure. |
|
Parses CDIP data from a web request. |
|
Returns the name of the predefined region in which the given coordinates reside. |
|
Returns data from the WPTO wave hindcast hosted on AWS at the specified latitude and longitude point(s). |
- mhkit.wave.IO.CDIP.cdip_request_parse_workflow(options)
CDIP_REQUEST_PARSE_WORKFLOW Parses CDIP data from a web request
Requests data for a station number (from http://cdip.ucsd.edu/) and parses. Years may be non-consecutive e.g. [2001, 2010]. Time may be sliced by dates (start_date or end date in YYYY-MM-DD). By default 2D variables are not parsed if all 2D variables are needed.
- Parameters:
station_number (
string) – Station number of CDIP wave buoyparameters (
string or array of strings) – Parameters to return. If nan, will return all variables except 2D-variables.years (
int or array of int) – Year date, e.g. 2001 or [2001, 2010]start_date (
string) – Start date in YYYY-MM-DD, e.g. ‘2012-04-01’end_date (
string) – End date in YYYY-MM-DD, e.g. ‘2012-04-30’data_type (
string) – Either ‘historic’ or ‘realtime’all_2D_variables (
boolean) – Will return all 2D data. Enabling this will add significant processing time. If all 2D variables are not needed it is recommended to pass 2D parameters of interest using the ‘parameters’ keyword and leave this set to the default false.
- Returns:
data (structure) –
- ‘data’: structure
grouped 1D and 2D structures of array data with datetimes
- ’metadata’: structure
Anything not of length time, including buoy name
- mhkit.wave.IO.hindcast.request_wpto(data_type, parameter, lat_lon, year, api_key)
Returns data from the WPTO wave hindcast hosted on AWS at the specified latitude and longitude point(s), or the closest available pont(s). Visit https://registry.opendata.aws/wpto-pds-us-wave/ for more information about the dataset and available locations and years. NOTE: To access the WPTO hindcast data, you will need to configure h5pyd for data access on HSDS. Please see the WPTO_hindcast_example notebook for more information.
- Parameters:
data_type (
string) – Data set type of interest Options: “3-hour” “1-hour”parameter (
string or list of strings) –Dataset parameter to be downloaded 3-hour dataset options: ‘directionality_coefficient’, ‘energy_period’, ‘maximum_energy_direction’
’mean_absolute_period’, ‘mean_zero-crossing_period’, ‘omni-directional_wave_power’, ‘peak_period’ ‘significant_wave_height’, ‘spectral_width’, ‘water_depth’
- 1-hour dataset options: ‘directionality_coefficient’, ‘energy_period’, ‘maximum_energy_direction’
’mean_absolute_period’, ‘mean_zero-crossing_period’, ‘omni-directional_wave_power’, ‘peak_period’ ‘significant_wave_height’, ‘spectral_width’, ‘water_depth’, ‘maximim_energy_direction’, ‘mean_wave_direction’, ‘frequency_bin_edges’, ‘directional_wave_spectrum’
lat_lon (
tuple or list of tuples) – Latitude longitude pairs at which to extract data.year (
float) – Year to be accessed. The years 1979-2010 available.api_key – API key obtained from https://developer.nrel.gov/signup/
- mhkit.wave.IO.hindcast.region_selection(lat_lon, data_type)
Returns the name of the predefined region in which the given coordinates reside. Can be used to check if the passed lat/lon pair is within the WPTO hindcast dataset.
- Parameters:
lat_lon (
matrix) – Latitude (first column) and longitude (second column) coordinates as numericsdata_type (
string) – Data set type of interest Options: “3-hour” “1-hour”
- Returns:
region (string) – Name of predefined region for given coordinates
- mhkit.wave.IO.NDBC.read_NDBC_file(file_name, varargin)
Reads a NDBC wave buoy data file (from https://www.ndbc.noaa.gov) into a structure.
Realtime and historical data files can be loaded with this function.
Note: With realtime data, missing data is denoted by “MM”. With historical data, missing data is denoted using a variable number of # 9’s, depending on the data type (for example: 9999.0 999.0 99.0). ‘N/A’ is automatically converted to missing data.
Data values are converted to float/int when possible. Column names are also converted to float/int when possible (this is useful when column names are frequency).
- Parameters:
file_name (
string) – Name of NDBC wave buoy data filemissing_value (
vector of values (optional)) – vector of values that denote missing data
- Returns:
data (Structure)
data.Data: named according to header row
data.time: given in datetime
data.units: the units for each data entry
OR if a spectra data NDBC file
data.spectrum: spectra data
data.time: given in datetime
data.frequency: spectral frequency
- mhkit.wave.IO.NDBC.NDBC_available_data(parameter, options)
For a given parameter this will return a structure of years, station IDs and file names that contain that parameter data.
- Parameters:
parameter (
string) – ‘swden’ : ‘Raw Spectral Wave Current Year Historical Data’ ‘stdmet’: ‘Standard Meteorological Current Year Historical Data’buoy_number (
string (optional)) – Buoy Number. 5-character alpha-numeric station identifier to call: NDBC_available_data(parameter,”buoy_number”,buoy_number)proxy (
None) –Parameter is now deprecated. To request data from behind a firewall, configure in MATLAB Preferences by navigating to:
Home -> Environment -> Preferences
- then:
MATLAB -> Web -> Use a proxy server to connect to the Internet
- See the following for details:
https://www.mathworks.com/help/matlab/import_export/proxy.html
- Returns:
available_data (structure) – Structure with station_id, years, and NDBC file names.
- mhkit.wave.IO.NDBC.NDBC_request_data(parameter, filenames)
Requests data by filenames and returns a structure of structures for each filename passed. If filenames for a single buoy are passed then the yearly structures in the returned structure (ndbc_data) are indexed by year (e.g. ndbc_data.year_2004). If multiple buoy ids are passed then the returned dictionary is indexed by buoy id and year (e.g. ndbc_data[‘46022’][‘2014’]).
- Parameters:
parameter (
string) – ‘swden’ : ‘Raw Spectral Wave Current Year Historical Data’ ‘stdmet’ : ‘Standard Meteorological Current Year Historical Data’filenames (
array of strings) – Data filenames on https://www.ndbc.noaa.gov/data/historical/{parameter}/proxy (
None) –Parameter is now deprecated. To request data from behind a firewall, configure in MATLAB Preferences by navigating to:
Home -> Environment -> Preferences
- then:
MATLAB -> Web -> Use a proxy server to connect to the Internet
- See the following for details:
https://www.mathworks.com/help/matlab/import_export/proxy.html
- Returns:
ndbc_data (Structure) – Structure of structures broken down by buoy and years of data.
- mhkit.wave.IO.SWAN.swan_read_block(swan_file)
Reads in SWAN ASCII block format output and returns a data structure containing modeled values and assocuated metadata
- Parameters:
swan_file (
string) – SWAN file name to import- Returns:
data (Structure)
- mhkit.wave.IO.SWAN.swan_read_table(swan_file)
Reads in SWAN ASCII table format output and returns a Matlab structure with modeled data and assocaited metadata.
- Parameters:
swan_file (
string) – SWAN file name to import- Returns:
data (Structure)
Resource
The resource submodule contains methods to compute wave energy spectra and various metrics from the spectra.
The following options exist to compute wave energy spectra:
Functions |
Description |
|---|---|
|
Calculates Jonswap spectrum from wave data. |
|
Calculates Pierson Moskowitz spectrum from wave data. |
|
Calculates wave spectra from wave probe timeseries. |
The following metrics can be computed from the spectra:
Functions |
Description |
|---|---|
|
Calculate the average creat period from spectra. |
|
Calculates the average wave period from spectra |
|
Calculates wave average zero crossing period from spectra |
|
Calculates depth regime based on wavelength and height |
|
Calculates the omnidirectional wave energy flux of the spectra |
|
Calculates the energy period |
|
Calculates environmental contours of extreme sea states |
|
Calculates the Nth frequency moment of the spectrum |
|
Calculates wave energy period from spectra |
|
Get Hs points along a specified environmental contour using user-defined T values. |
|
Sample a sea state between contours of specified return periods. |
|
Calculates wave height from spectra |
|
Calculates bandwidth from spectra |
|
Calculates wave spectral width from spectra |
|
Calculates wave elevation time series from spectrum using a random phase |
|
Calculates wave celerity (group velocity) |
|
Calculates wave length from wave number |
|
Calculates wave number |
- mhkit.wave.resource.spectral_bandwidth(S, varargin)
Calculates bandwidth from spectra
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- frequency_bins: vector (optional)
Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
e (double) – Spectral BandWidth
- mhkit.wave.resource.samples_contour(t_samples, t_contour, hs_contour)
Get Hs points along a specified environmental contour using user-defined T values.
- Parameters:
t_samples (
array) – Points for sampling along return contourt_contour (
array) – T values along contourhs_contour (
array) – Hs values along contour
- Returns:
hs_samples (array) – points sampled along return contour
- mhkit.wave.resource.spectral_width(S, varargin)
Calculates wave spectral width from spectra
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- frequency_bins: vector (optional)
Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
e0 (float) – Spectral Width
- mhkit.wave.resource.significant_wave_height(S, varargin)
Calculates significant wave height Hm0 from spectra
- Parameters:
S (
structure or numeric array) –- If structure:
S.spectrum: Spectral density [m^2/Hz] S.frequency: frequency [Hz]
- If numeric:
S is assumed to be spectral density vector/matrix varargin{1} must contain frequency vector
frequency_bins (
vector (optional)) – Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
Hm0 (double) – Significant wave height [m] index by S.columns
- mhkit.wave.resource.elevation_spectrum(ts, sample_rate, nnft, time, options)
Calculates wave spectra from wave probe timeseries
- Parameters:
ts (
matrix or table) – Wave probe time-series data, with each column a different time seriessample_rate (
float) – Data frequency (Hz)nnft (
integer) – Number of bins in the Fast Fourier Transformtime (
vector or table) – time (s)window (
string (Optional)) – Signal window type. “hamming” is used by default given the broadband nature of waves. See scipy.signal.get_window for more options. to call: elevation_spectrum(ts,sample_rate,nnft,time,”window”,window)detrend (
logical (Optional)) – Specifies if a linear trend is removed from the data before calculating the wave energy spectrum. Data is detrended by default. to call: elevation_spectrum(ts,sample_rate,nnft,time,”detrend”,detrend)noverlap (
integer (Optional)) – Number of points to overlap between segments. If None, noverlap = nperseg / 2. Defaults to None. to call: elevation_spectrum(ts,sample_rate,nnft,time,”noverlap”,noverlap)
- Returns:
S (structure)
S.spectrum: vector or matrix Spectral Density (m^2/Hz) per probe
S.type: ‘Spectra from Time Series’
S.frequency: frequency [Hz]
S.sample_rate: sample_rate
S.nnft: nnft
- mhkit.wave.resource.average_crest_period(S, varargin)
Calculates the average crest period
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- frequency_bins: vector (optional)
Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
Tavg (double) – Average crest Period (s)
- mhkit.wave.resource.energy_flux(S, h, options)
- Parameters:
S (
Spectral Density (m^2/Hz)) –- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
h (
double) – Water depth (m)deep (
logical (optional)) – If True use the deep water approximation. Default False. When False a depth check is run to check for shallow water. The ratio of the shallow water regime can be changed using the ratio keyword. to call: energy_flux(S,h,”deep”,py.True)rho (
double (optional)) – water density (kg/m^3) to call: energy_flux(S,h,”rho”,rho)g (
double (optional)) – gravitational acceleration (m/s^2) to call: energy_flux(S,h,”g”,g)ratio (
double or int (optional)) – Only applied if depth=False. If h/l > ratio, water depth will be set to deep. Default ratio = 2. to call: energy_flux(S,h,”ratio”,1.5)
- Returns:
J (double) – Omni-directional wave energy flux (W/m)
- mhkit.wave.resource.samples_full_seastate(x1, x2, points_per_interval, return_periods, sea_state_duration, method, bin_size)
Sample a sea state between contours of specified return periods.
This function is used for the full sea state approach for the extreme load. See Coe et al. 2018 for more details. It was originally part of WDRT.
Coe, R. G., Michelen, C., Eckert-Gallup, A., & Sallaberry, C. (2018). Full long-term design response analysis of a wave energy converter. Renewable Energy, 116, 356-366.
- x1array
Component 1 data.
- x2array
Component 2 data.
- point_per_intervalint
Number of sample points to be calculated per contour interval.
- return_periodsarray
Vector of return periods that define the contour intervals in which samples will be taken. Values must be greater than zero and must be in increasing order.
- sea_state_durationdouble
x1 and x2 sample rate (seconds)
- methodstr or list
Copula method to apply. Currently only ‘PCA’ is implemented.
- bin_sizeint
Number of data points in each bin (250 is recommended).
- Hs_Samplesarray
Vector of Hs values for each sample point.
- Te_Samplesarray
Vector of Te values for each sample point.
- weight_pointsarray
Vector of probabilistic weights for each sampling point to be used in risk calculations.
- mhkit.wave.resource.depth_regime(l, h, options)
Calculates the depth regime based on wavelength and height Deep water: h/l > ratio This function exists so sinh in wave celerity doesn’t blow up to infinity.
P.K. Kundu, I.M. Cohen (2000) suggest h/l >> 1 for deep water (pg 209) Same citation as above, they also suggest for 3% accuracy, h/l > 0.28 (pg 210) However, since this function allows multiple wavelengths, higher ratio numbers are more accurate across varying wavelengths.
- Parameters:
l (
vector) – wave length (m)h (
double) – Water depth (m)ratio (
double or int (optional)) – If h/l > ratio, water depth will be set to deep. Default ratio = 2 to call: energy_flux(k,h,”ratio”,ratio)
- Returns:
depth_reg (boolean or boolean array) – Boolean True if deep water, False otherwise
- mhkit.wave.resource.wave_celerity(k, h, options)
Calculates wave celerity (group velocity)
- Parameters:
k (
wave number (1/m)) –- structure of form:
k.values= wave number
k.frequency= frequency (Hz)
h (
double) – Water depth (m)g (
double (optional)) – gravitational acceleration (m/s^2) to call: energy_flux(k,h,”g”,g)depth_check (
bool (optional)) – If True check depth regime. Default False. to call: energy_flux(k,h,”depth_check”,py.True)ratio (
double or int (optional)) – Only applied if depth_check=True. If h/l > ratio, water depth will be set to deep. Default ratio = 2 to call: energy_flux(k,h,”ratio”,ratio)
- Returns:
Cg (structure)
Cg.values: water celerity
Cg.frequency [Hz]
Cg.h: height [m]
- mhkit.wave.resource.wave_number(f, h, options)
Calculates wave number
- Parameters:
f (
frequency (Hz)) – vector or numpy arrayh (
float) – Water depth (m)rho (
float (optional)) – water density (kg/m^3) to call: wave_number(f,h,”rho”,rho)g (
float (optional)) – gravitational acceleration (m/s^2) to call: wave_number(f,h,”g”,g)
- Returns:
k (structure)
k.values: wave number
k.frequency: frequency [Hz]
- mhkit.wave.resource.wave_length(k)
Calculates wave length from wave number To compute: 2*pi/wavenumber
- Parameters:
k (
wave number (1/m)) – intiger, double, or vector- Returns:
l (double or array) – Wave length [m] indexed by frequency
- mhkit.wave.resource.pierson_moskowitz_spectrum(frequency, Tp, Hs)
Calculates Pierson-Moskowitz Spectrum from IEC TS 62600-2 ED2 Annex C.2 (2019)
- Parameters:
frequency (
vector) – Wave frequency [Hz]Tp (
float) – Peak period [s]Hs (
float) – Significant wave height [m]
- Returns:
S (structure with fields) – .spectrum = Spectral density [m^2/Hz] .frequency = Frequency [Hz] .type = ‘Pierson-Moskowitz (Hs,Tp)’
- mhkit.wave.resource.standardize_wave_spectra_frequency(frequency, spectrum, varargin)
Standardize frequency and spectrum data for wave calculations
This utility function standardizes frequency vectors, spectra, and frequency bins following MHKiT-Python conventions for consistent numerical integration.
- Parameters:
frequency (
vector) – Frequency vector (Hz)spectrum (
vector or matrix) – Spectral density datafreq_bins (
vector (optional)) – Bin widths for frequency. If not provided, calculated from frequency differences
- Returns:
frequency (vector) – Filtered frequency vector (> 1e-12 Hz), column format
spectrum (vector or matrix) – Filtered spectrum data matching frequency
freq_bins (vector) – Frequency bin widths, column format
- mhkit.wave.resource.environmental_contours(x1, x2, dt, period, method, options)
Calculates environmental contours of extreme sea states using the improved joint probability distributions with the inverse first-order reliability method (IFORM) probability for the desired return period (period). Given the period of interest a circle of iso-probability is created in the in the PCA joint probability (x1, x2) reference frame. Using the joint probability value the CDF of the marginal distribution is used to find the quantile of each component. Finally, using the improved PCA methodology the component 2 contour lines are calculated from component 1 using the relationships defined in Exkert-Gallup et. al. 2016.
Eckert-Gallup, A. C., Sallaberry, C. J., Dallman, A. R., & Neary, V. S. (2016). Application of principal component analysis (PCA) and improved joint probability distributions to the inverse first-order reliability method (I-FORM) for predicting extreme sea states. Ocean Engineering, 112, 307-319.
- Parameters:
x1 (
vector) – component 1 datax2 (
vector) – component 2 datadt (
double) – x1 and x2 sample rate (seconds)period (
scalar or vector) – Desired return period (years) for calculation of environmental contour, can be a scalar or a vector.method (
string or list) – Copula method to apply. Options include [‘PCA’,’gaussian’, ‘gumbel’, ‘clayton’, ‘rosenblatt’, ‘nonparametric_gaussian’, ‘nonparametric_clayton’, ‘nonparametric_gumbel’, ‘bivariate_KDE’ ‘bivariate_KDE_log’]**OPTIONS –
- min_bin_count: int
Passed to _copula_parameters to sets the minimum number of bins allowed. Default = 40.
- initial_bin_max_val: int, double
Passed to _copula_parameters to set the max value of the first bin. Default = 1.
- bin_val_size: int, double
Passed to _copula_parameters to set the size of each bin after the initial bin. Default 0.25.
- nb_steps: int
Discretization of the circle in the normal space is used for copula component calculation. Default nb_steps=1000.
- bandwidth:
Must specify bandwidth for bivariate KDE method. Default = None.
- Ndata_bivariate_KDE: int
Must specify bivariate KDE method. Defines the contoured space from which samples are taken. Default = 100.
- max_x1: double
Defines the max value of x1 to discretize the KDE space
- max_x2: double
Defines the max value of x2 to discretize the KDE space
- PCA: dict
If provided, the principal component analysis (PCA) on x1, x2 is skipped. The PCA will be the same for a given x1, x2 therefore this step may be skipped if multiple calls to environmental contours are made for the same x1, x2 pair. The PCA dict may be obtained by setting return_fit=True when calling the PCA method.
- return_fit: boolean
Will return fitting parameters used for each method passed. Default False.
- Returns:
environmental_contour (Structure) – Structure with fields contour1, contour2, and optionally PCA
- mhkit.wave.resource.average_zero_crossing_period(S, varargin)
Calculates wave average zero crossing period from spectra
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra,x)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- frequency_bins: vector (optional)
Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
Tz (double) – Average Zero Crossing Period (s)
- mhkit.wave.resource.peak_period(S)
Calculates wave energy period from spectra
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- Returns:
Tp float – Wave Peak Period (s)
- mhkit.wave.resource.average_wave_period(S, varargin)
Calculates the average wave period
- Parameters:
S (
Spectral Density (m^2/Hz)) –- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(frequency,spectra)
OR
- structure of form:
S.spectrum: Spectral Density (m^2/Hz)
S.type: String of the spectra type, i.e. Bretschneider, time series, date stamp etc.
S.frequency: frequency (Hz)
- frequency_bins: vector (optional)
Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
Tavg (float) – Mean wave period (s)
- mhkit.wave.resource.surface_elevation(S, time_index, options)
Calculates wave elevation time-series from spectrum
- Parameters:
S (
structure) –- Spectral data with fields:
S.spectrum: (n_freq x 1) spectral density [m^2/Hz] S.frequency: (n_freq x 1) frequency vector [Hz]
time_index (
vector) – Time used to create the wave elevation time-series [s], for example, time = 0:0.01:100options (
structure (optional)) –- Optional fields:
seed: Random seed (default = 123) frequency_bins: Bin widths for frequency of S. Required for unevenly sized bins phases: Explicit phases for frequency components (overrides seed)
for example, phases = rand(length(S.frequency),1) * 2 * pi
- method: Method used to calculate the surface elevation. ‘ifft’
(Inverse Fast Fourier Transform) used by default if the given frequency_bins==[] or is evenly spaced. ‘sum_of_sines’ explicitly sums each frequency component and used by default if uneven frequency_bins are provided. The ‘ifft’ method is significantly faster.
- Returns:
wave_elevation (structure) –
- Generated wave elevation with fields:
elevation: Wave surface elevation [m] time: time vector [s] type: description string
- mhkit.wave.resource.jonswap_spectrum(frequency, Tp, Hs, gamma)
Calculates JONSWAP Spectrum from IEC TS 62600-2 ED2 Annex C.2 (2019)
- Parameters:
frequency (
vector) – Frequency [Hz]Tp (
float) – Peak period [s]Hs (
float) – Significant wave height [m]gamma (
float (optional)) – Peak enhancement factor
- Returns:
S (structure with fields) – .spectrum = Spectral density [m^2/Hz] .frequency = Frequency [Hz] .type = ‘JONSWAP (Hs,Tp)’
- mhkit.wave.resource.energy_period(S, varargin)
Calculate wave energy period (Te) seconds from power spectral density (PSD)
- Parameters:
S (
structure or numeric array) –- If structure:
S.spectrum - Spectral density [m^2/Hz] S.frequency - Frequency [Hz]
- If numeric:
S is spectral density array (vector or matrix) varargin{1} must contain frequency vector
frequency_bins (
vector (optional)) – Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
Te (double) – Wave energy period [seconds] indexed by S.columns
- mhkit.wave.resource.frequency_moment(S, N, varargin)
Calculates the Nth frequency moment of the spectrum
- Parameters:
S (
structure or numeric array) –- If structure:
S.spectrum - Spectral density [m^2/Hz] S.frequency - Frequency [Hz]
- If numeric:
S is spectral density array (vector or matrix) varargin{1} must contain frequency vector
N (
int) – Moment (0 for 0th, 1 for 1st ….)frequency_bins (
vector (optional)) – Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
m (double) – Nth Frequency Moment indexed by S.columns
Performance
The performance submodule contains functions to compute capture length, statistics, performance matrices, and mean annual energy production.
Functions |
Description |
|---|---|
|
Calculates the capture length (often called capture width). |
|
Generates a capture length matrix for a given statistic |
|
Calculates mean annual energy production (MAEP) from matrix data along with data frequency in each bin |
|
Calculates mean annual energy production (MAEP) from timeseries |
|
Generates a power matrix from a capture length matrix and wave energy flux matrix |
|
Generates a wave eneergy flux matrix for a given statistic |
|
High-level function to compute power performance quantities of interest following IEC TS 62600-100 for given wave spectra. |
- mhkit.wave.performance.capture_length(Power, J)
Calculates the capture length (often called capture width).
- Parameters:
P (
array or vector) – Power [W]J (
array or vector) – Omnidirectional wave energy flux [W/m]
- Returns:
L (vector) – Capture length [m]
- mhkit.wave.performance.wave_energy_flux_matrix(Hm0, Te, J, statistic, Hm0_bins, Te_bins)
Generates a wave eneergy flux matrix for a given statistic
Note that IEC/TS 62600-100 requires capture length matrices for the mean, std, count, min, and max.
- Parameters:
Hm0 (
numpy array or vector) – Significant wave height from spectra [m]Te (
numpy array or vector) – Energy period from spectra [s]J (
numpy array or vector) – wave energy flux from spectra [W/m]statistic (
string) – Statistic for each bin, options include: ‘mean’, ‘std’, ‘median’, ‘count’, ‘sum’, ‘min’, ‘max’, and ‘frequency’. Note that ‘std’ uses a degree of freedom of 1 in accordance with IEC/TS 62600-100.Hm0_bins (
numpy array or vector) – Bin centers for Hm0 [m]Te_bins (
numpy array or vector) – Bin centers for Te [s]
- Returns:
WEFM (Structure)
WEFM.values
WEFM.stat
WEFM.Hm0_bins
WEFM.Te_bins
- mhkit.wave.performance.mean_annual_energy_production_matrix(LM, JM, frequency)
Calculates mean annual energy production (MAEP) from matrix data along with data frequency in each bin
- Parameters:
LM –
- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(Hm0_bins,L)
OR
structure of form:
LM.values
LM.stat
LM.Hm0_bins
LM.Te_bins
- Returns:
maep (float) – Mean annual energy production
- mhkit.wave.performance.power_matrix(LM, JM)
Generates a power matrix from a capture length matrix and wave energy flux matrix
- Parameters:
LM –
- Pandas data frame
To make a pandas data frame from user supplied frequency and spectra use py.mhkit_python_utils.pandas_dataframe.spectra_to_pandas(Hm0_bins,L)
OR
structure of form:
LM.values
LM.stat
LM.Hm0_bins
LM.Te_bins
- Returns:
PM (Structure)
PM.values: Power matrix
PM.stat: statistic of the matrix (i.e. “mean”, “max”, etc.) (string)
PM.Hm0_bins
PM.Te_bins
- mhkit.wave.performance.capture_length_matrix(Hm0, Te, L, statistic, Hm0_bins, Te_bins)
Generates a capture length matrix for a given statistic
Note that IEC/TS 62600-100 requires capture length matrices for the mean, std, count, min, and max.
- Parameters:
Hm0 (
numpy array or vector) – Significant wave height from spectra [m]Te (
numpy array or vector) – Energy period from spectra [s]L (
numpy array or vector) – Capture length [m]statistic (
string) – Statistic for each bin, options include: ‘mean’, ‘std’, ‘median’, ‘count’, ‘sum’, ‘min’, ‘max’, and ‘frequency’. Note that ‘std’ uses a degree of freedom of 1 in accordance with IEC/TS 62600-100. or a callable function of python typeHm0_bins (
numpy array or vector) – Bin centers for Hm0 [m]Te_bins (
numpy array or vector) – Bin centers for Te [s]
- Returns:
clm (structure)
clm.values
clm.stat
clm.Hm0_bins
clm.Te_bins
- mhkit.wave.performance.mean_annual_energy_production_timeseries(L, J)
Calculates mean annual energy production (MAEP) from timeseries
- Parameters:
L (
numpy array or vector) – Capture lengthJ (
numpy array or vector) – Wave energy flux
- Returns:
maep (float) – Mean annual energy production
- mhkit.wave.performance.power_performance_workflow(S, h, P, statistic, options)
High-level function to compute power performance quantities of interest following IEC TS 62600-100 for given wave spectra.
- Parameters:
S (
structure with fields:) – S.spectrum: Spectral Density [m^2/Hz] S.frequency: frequency [Hz] S.time : time [datetime]h (
integer) – Water depth [m]P (
array or vector) – Power [W]statistic (
string or array of strings) – Capture length statistics for plotting options include: “mean”, “std”, “median”, “count”, “sum”, “min”, “max”, and “frequency”. Note that “std” uses a degree of freedom of 1 in accordance with IEC/TS 62600-100. To output capture length matrices for multiple binning parameters, define as a string array: statistic = [“”, “”, “”];savepath (
string (optional)) – Path to save figure. to call: power_performance_wave(S,h,P,statistic,”savepath”,savepath)rho (
float (optional)) – Water density [kg/m^3] to call: power_performance_wave(S,h,P,statistic,”rho”,rho)g (
float (optional)) – Gravitational acceleration [m/s^2] to call: power_performance_wave(S,h,P,statistic,”g”,g)frequency_bins (
vector (optional)) – Bin widths for frequency of S. Required for unevenly sized bins
- Returns:
cl_matrix (figure) – Capture length matrix
maep_matrix (float) – Mean annual energy production
Graphics
The :graphics submodule contains functions to plot wave data and related metrics.
Functions |
Description |
|---|---|
|
Plots wave elevation timeseries |
|
Plots an overlay of the x1 and x2 variables to the calculated environmental contours. |
|
Plots the matrix with Hm0 and Te on the y and x axis |
|
Plots wave amplitude spectrum |
|
Plots, in the style of Chakrabarti (2005), relative importance of viscous,inertia, and diffraction phemonena |
|
Creates monthl average boxplots of significant wave height |
|
Creates subplots of environmental resource from cdip data |
- mhkit.wave.graphics.plot_compendium(Hs, Tp, Dp, time, options)
PLOT_COMPENDIUM Creates subplots of wave height, peak period and direction Creates subplots showing: Significant Wave Height (Hs), Peak Period (Tp), and Direction (Dp). Developed based on: http://cdip.ucsd.edu/themes/media/docs/documents/html_pages/compendium.html
- Parameters:
Hs (
double array) – significant wave heightTp (
double array) – peak periodDp (
double array) – directiontime (
datetime array) – timestampsbuoy_title (
string (optional)) – figure suptitle (super title)
- Returns:
f (figure object)
- mhkit.wave.graphics.plot_boxplot(Hs, time, options)
PLOT_BOXPLOT Creates monthly-averaged boxplots of significant wave height Creates monthly-averaged boxplots of significant wave height (Hs) Developed based on: http://cdip.ucsd.edu/themes/media/docs/documents/html_pages/annualHs_plot.html
- Parameters:
Hs (
double array) – significant wave heighttime (
datetime array) – timestampsbuoy_title (
string (optional)) – figure suptitle (super title)
- Returns:
f (figure object)
- mhkit.wave.graphics.plot_matrix(M, Mtype, options)
Plots the matrix with Hm0 and Te on the y and x axis
- Parameters:
M (
structure) –M.values: matrix
M.Hm0_bins
M.Te_bins
M.stat
Mtype (
string) –- type of matrix (i.e. power, capture length, etc.) to be used
in plot title
- options: name-value pairs
- savepath: string (optional)
path and filename to save figure.
- annotate: logical (optional)
toggle text annotations on/off (default: true)
to call: plot_matrix(M, Mtype,”savepath”,savepath,”annotate”,false)
- Returns:
figure (plot of the matrix)
- mhkit.wave.graphics.plot_chakrabarti(H, lambda_w, D, options)
Plots, in the style of Chakrabarti (2005), relative importance of viscous, inertia, and diffraction phemonena
Chakrabarti, Subrata. Handbook of Offshore Engineering (2-volume set). Elsevier, 2005.
- Parameters:
H (
integer, double or vector) – Wave height [m]lambda_w (
integer, double or vector) – Wave length [m]D (
integer, double or vector of) – Characteristic length [m]savepath (
string (optional)) – path and filename to save figure. to call: plot_chakrabarti(H,lambda_w,D,”savepath”,savepath)
- Returns:
figure (figure)
Plots wave force regime as Keulegan-Carpenter parameter versus
diffraction parameter
- Examples
Using Integers >> D = 5 >> H = 8 >> lambda_w = 200 >> plot_chakrabarti(H,lambda_w,D)
Using vector >> D = linspace(5,15,5) >> H = 8*ones(size(D)) >> lambda_w = 200*ones(size(D)) >> plot_chakrabarti(H,lambda_w,D)
- mhkit.wave.graphics.plot_environmental_contours(x1, x2, x1_contour, x2_contour, options)
Plots an overlay of the x1 and x2 variables to the calculated environmental contours.
- Parameters:
x1 (
vector) – x-axis datax2 (
vector) – y-axis datax1_countour (
Table) – Calculated x1 contour valuesx2_countour (
Table) – Calculated x2 contour valuesx_label (
string (optional)) – x-axis label. Default None. to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”x_label”,x_label)y_label (
string (optional)) – y-axis label. Default None. to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”y_label”,y_label)data_label (
string (optional)) – Legend label for x1, x2 data (e.g. ‘Buoy 46022’). Default None. to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”data_label”,data_label)contour_label (
string or array of strings (optional)) – Legend label for x1_contour, x2_contour countor data (e.g. ‘100-year contour’). Default None. to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”contour_label”,contour_label)title (
string (optional)) – title for the plot to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”title”,title)savepath (
string (optional)) – path and filename to save figure. to call: plot_environmental_contours(x1,x2,x1_contour,x2_contour,”savepath”,savepath)
- Returns:
figure (figure) – Envoronmental contour plot
- mhkit.wave.graphics.plot_elevation_timeseries(wave_elevation, options)
Plots wave elevation timeseries
- Parameters:
wave_elevation (
Structure of the following form:) –wave_elevation.elevation: elevation [m]
wave_elevation.time: time (s);
title (
string (optional)) – title for the plot to call: plot_elevation_timeseries(wave_elevation,”title”,title)savepath (
string (optional)) – path and filename to save figure. to call: plot_elevation_timeseries(wave_elevation,”savepath”,savepath)
- Returns:
figure (figure) – Plot of wave elevation vs. time
- mhkit.wave.graphics.plot_spectrum(wave_spectra, options)
Plots wave amplitude spectrum
- Parameters:
wave_spectra
- Returns:
figure (figure)
Plot of wave amplitude spectra versus omega