bluemath_tk.waves package
Submodules
bluemath_tk.waves.binwaves module
bluemath_tk.waves.calibration module
bluemath_tk.waves.climate module
bluemath_tk.waves.estela module
bluemath_tk.waves.greenswell module
bluemath_tk.waves.partitioning module
bluemath_tk.waves.series module
- bluemath_tk.waves.series.energy_spectrum(hs: float, tp: float, gamma: float, duration: int) List[float] [source]
Calculate JONSWAP energy spectrum for given wave parameters.
This function computes the energy spectrum S(f) using the JONSWAP spectral model, which is commonly used for wind-generated waves in deep water.
- Parameters:
hs (float) – Significant wave height in meters.
tp (float) – Peak period in seconds.
gamma (float) – JONSWAP peak enhancement factor (typically 1-7, default 3.3).
duration (int) – Number of frequency points for the spectrum.
- Returns:
Energy spectrum values S(f) for each frequency.
- Return type:
List[float]
Notes
The JONSWAP spectrum combines the Pierson-Moskowitz spectrum with a peak enhancement factor. The spectrum is defined as:
S(f) = Beta * (Hs^2) * (Tp^-4) * (f^-5) * exp(-1.25 * (Tp*f)^-4) * gamma^exp(-((Tp*f-1)^2)/(2*sigma^2))
where sigma varies based on frequency relative to peak frequency.
- bluemath_tk.waves.series.series_Jonswap(waves: Dict[str, float | int]) ndarray [source]
Generate surface elevation time series from JONSWAP spectrum.
This function generates a time series of surface elevation using the JONSWAP spectral model. The elevation is computed by summing harmonic components with random phases.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters: - ‘H’: Significant wave height (m) - ‘T’: Peak period (s) - ‘gamma’: JONSWAP peak enhancement factor - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘tendc’: Simulation period (s)
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The surface elevation is generated by: 1. Computing the JONSWAP energy spectrum 2. Converting spectrum to amplitude components 3. Summing harmonic components with random phases
- bluemath_tk.waves.series.series_Jonswap_bimodal(waves: Dict[str, float | int]) ndarray [source]
Generate surface elevation time series from bimodal JONSWAP spectrum.
This function generates a time series of surface elevation using a bimodal JONSWAP spectral model, combining two wave systems with different characteristics.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters for two wave systems: - ‘Hs1’: Significant wave height of system 1 (m) - ‘Tp1’: Peak period of system 1 (s) - ‘gamma1’: JONSWAP peak enhancement factor for system 1 - ‘Hs2’: Significant wave height of system 2 (m) - ‘Tp2’: Peak period of system 2 (s) - ‘gamma2’: JONSWAP peak enhancement factor for system 2 - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘tendc’: Simulation period (s)
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The bimodal spectrum is created by summing two JONSWAP spectra, each representing a different wave system (e.g., wind sea and swell).
- bluemath_tk.waves.series.series_TMA(waves: Dict[str, float | int], depth: float) ndarray [source]
Generate surface elevation time series from TMA spectrum.
This function generates a time series of surface elevation using the TMA spectral model, which extends the JONSWAP spectrum to include depth effects.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters: - ‘H’: Significant wave height (m) - ‘T’: Peak period (s) - ‘gamma’: JONSWAP peak enhancement factor - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘comptime’: Simulation period (s)
depth (float) – Water depth in meters.
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The TMA spectrum modifies the JONSWAP spectrum to account for depth effects using a depth function fi that depends on wave number and depth.
- bluemath_tk.waves.series.series_TMA_bimodal(waves: Dict[str, float | int], depth: float) ndarray [source]
Generate surface elevation time series from bimodal TMA spectrum.
This function generates a time series of surface elevation using a bimodal TMA spectral model, combining two wave systems with different characteristics.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters for two wave systems: - ‘Hs1’: Significant wave height of system 1 (m) - ‘Tp1’: Peak period of system 1 (s) - ‘gamma1’: JONSWAP peak enhancement factor for system 1 - ‘Hs2’: Significant wave height of system 2 (m) - ‘Tp2’: Peak period of system 2 (s) - ‘gamma2’: JONSWAP peak enhancement factor for system 2 - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘tendc’: Simulation period (s)
depth (float) – Water depth in meters.
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The bimodal spectrum is created by summing two TMA spectra, each representing a different wave system (e.g., wind sea and swell).
- bluemath_tk.waves.series.series_regular_bichromatic(waves: Dict[str, float | int]) ndarray [source]
Generate bichromatic regular wave time series.
This function generates a time series of surface elevation for bichromatic regular waves, combining two monochromatic components.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters: - ‘T1’: Period of component 1 (s) - ‘T2’: Period of component 2 (s) - ‘H’: Wave height (m) - ‘WL’: Water level (m) - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘tendc’: Simulation period (s)
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The surface elevation is the sum of two cosine components: η(t) = (H/2) * cos(2πt/T1) + (H/2) * cos(2πt/T2)
- bluemath_tk.waves.series.series_regular_monochromatic(waves: Dict[str, float | int]) ndarray [source]
Generate monochromatic regular wave time series.
This function generates a time series of surface elevation for monochromatic regular waves with constant amplitude and frequency.
- Parameters:
waves (Dict[str, Union[float, int]]) – Dictionary containing wave parameters: - ‘T’: Wave period (s) - ‘H’: Wave height (m) - ‘warmup’: Spin-up time (s) - ‘deltat’: Time step (s) - ‘comptime’: Simulation period (s)
- Returns:
2D array with shape (n_time_steps, 2) containing: - Column 0: Time values (s) - Column 1: Surface elevation values (m)
- Return type:
np.ndarray
Notes
The surface elevation follows a simple cosine function: η(t) = (H/2) * cos(2πt/T)
- bluemath_tk.waves.series.waves_dispersion(T: float, h: float) Tuple[float, float, float] [source]
Solve the wave dispersion relation to calculate wavelength, wave number, and phase speed.
This function solves the linear wave dispersion relation to determine wave properties in water of finite depth.
- Parameters:
T (float) – Wave period in seconds.
h (float) – Water depth in meters.
- Returns:
A tuple containing: - L: Wavelength in meters - k: Wave number in radians per meter - c: Phase speed in meters per second
- Return type:
Tuple[float, float, float]
Notes
The dispersion relation is solved iteratively using: L = (g*T²)/(2π) * tanh(2πh/L)
where g is gravitational acceleration (9.81 m/s²).
The wave number k = 2π/L and phase speed c = √(g*tanh(kh)/k).
bluemath_tk.waves.snakes module
bluemath_tk.waves.spectra module
- bluemath_tk.waves.spectra.spectral_analysis(water_level: ndarray, delta_time: float) Tuple[float, float, float, float] [source]
Performs spectral analysis of water level variable to separate wave components.
This function uses Welch’s method to estimate the power spectral density of the water level time series and then separates the energy into different frequency bands: incident waves (IC), infragravity waves (IG), and very low frequency waves (VLF). The significant wave heights are calculated for each component.
- Parameters:
water_level (np.ndarray) – Water level time series in meters. Should be a 1D array of measurements.
delta_time (float) – Time interval between consecutive samples in seconds.
- Returns:
A tuple containing the significant wave heights in meters: - Hsi: Total significant wave height - HIC: Incident wave significant height (f > 0.04 Hz) - HIG: Infragravity wave significant height (0.004 Hz < f < 0.04 Hz) - HVLF: Very low frequency wave significant height (0.001 Hz < f < 0.004 Hz)
- Return type:
Tuple[float, float, float, float]
Notes
The function uses the following frequency bands for wave component separation: - Incident waves (IC): f > 0.04 Hz - Infragravity waves (IG): 0.004 Hz < f < 0.04 Hz - Very low frequency waves (VLF): 0.001 Hz < f < 0.004 Hz
The significant wave height is calculated as H_s = 4 * sqrt(m0), where m0 is the zeroth moment of the spectrum.
Examples
>>> import numpy as np >>> water_level = np.random.randn(1000) # 1000 samples >>> delta_time = 1.0 # 1 second intervals >>> Hsi, HIC, HIG, HVLF = spectral_analysis(water_level, delta_time)
bluemath_tk.waves.statistics module
- bluemath_tk.waves.statistics.upcrossing(series: ndarray) Tuple[List[float], ndarray] [source]
Calculate wave heights and periods using the zero-crossing method.
This function performs wave analysis by detecting zero-crossings in a surface elevation time series. It removes the mean from the series and identifies upcrossings (where the signal crosses zero from negative to positive). For each wave cycle between upcrossings, it calculates the wave height (difference between maximum and minimum) and wave period (time between crossings).
- Parameters:
series (np.ndarray) – A 2D array with shape (n_samples, 2) where: - series[:, 0] contains the time values - series[:, 1] contains the surface elevation values
- Returns:
A tuple containing: - Hi : List of wave heights (amplitudes) for each wave cycle - Ti : Array of wave periods (time intervals between zero crossings)
- Return type:
Tuple[List[float], np.ndarray]
Notes
The function removes the mean from the surface elevation series before analysis
Wave heights are calculated as the difference between maximum and minimum values within each wave cycle
Wave periods are the time intervals between consecutive zero upcrossings
This method is commonly used in oceanography for wave analysis
Examples
>>> import numpy as np >>> time = np.linspace(0, 10, 1000) >>> elevation = np.sin(2 * np.pi * 0.5 * time) + 0.1 * np.random.randn(1000) >>> series = np.column_stack([time, elevation]) >>> heights, periods = upcrossing(series)
bluemath_tk.waves.superpoint module
- bluemath_tk.waves.superpoint.superpoint_calculation(stations_data: DataArray, stations_dimension_name: str, sectors_for_each_station: Dict[str, Tuple[float, float]]) DataArray [source]
Join multiple station spectral data for each directional sector.
- Parameters:
stations_data (xr.DataArray) – DataArray containing spectral data for multiple stations.
stations_dimension_name (str) – Name of the dimension representing different stations in the DataArray.
sectors_for_each_station (Dict[str, Tuple[float, float]]) – Dictionary mapping each station ID to a tuple of (min_direction, max_direction) representing the directional sector for that station.
- Returns:
A new DataArray where each point is the sum of spectral data from all stations for the specified directional sector.
- Return type:
xr.DataArray
Notes
If your stations_data is saved in different files, you can load them all and then concatenate them using xr.open_mfdataset function. Example below:
“path/to/station1.nc”, “path/to/station2.nc”, “path/to/station3.nc”
]
- def load_station_data(ds: xr.Dataset) -> xr.DataArray:
return ds.efth.expand_dims(“station”)
- stations_data = xr.open_mfdataset(
files, concat_dim=”station”, preprocess=load_station_data,
)
Module contents
Project: BlueMath_tk Sub-Module: waves Author: GeoOcean Research Group, Universidad de Cantabria Repository: https://github.com/GeoOcean/BlueMath_tk.git Status: Under development (Working)