proteus.WaveTools module

Tools for working with water waves.

The primary objective of this module is to provide solutions (exact and approximate) for the free surface deformation and subsurface velocity components of water waves. These can be used as boundary conditions, wave generation sources, and validation solutions for numerical wave codes.

class proteus.WaveTools.SteadyCurrent[source]

Bases: object

This class is used for generating a steady current

Parameters
  • U (numpy.ndarray) – Current velocity in vector form

  • mwl (float) – Still water level

  • rampTime (float) – Ramp time for current

eta(x, t)[source]

Calculates free surface elevation (SolitaryWave class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

mwl[source]
u(x, t)[source]

Calculates wave velocity vector (SolitaryWave class). :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

class proteus.WaveTools.SolitaryWave[source]

Bases: object

This class is used for generating 1st order solitary wave

Parameters
  • waveHeight (float) – Regular wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • g (numpy.ndarray) – Gravitational acceleration vector

  • waveDir (numpy.ndarray) – Wave direction in vector form

  • trans (numpy.ndarray) – Position vector of the peak

  • fast (bool) – Switch for optimised functions

c[source]
eta(x, t)[source]

Calculates free surface elevation (SolitaryWave class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

mwl[source]
u(x, t)[source]

Calculates wave velocity vector (SolitaryWave class). :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

class proteus.WaveTools.MonochromaticWaves[source]

Bases: object

This class is used for generating regular waves in both linear and nonlinear regimes. See Dean and Dalrymple 1994 for equations.

Parameters
  • period (float) – Regular wave period

  • waveHeight (float) – Regular wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • g (numpy.ndarray) – Gravitational acceleration vector

  • waveDir (numpy.ndarray) – Wave direction in vector form

  • wavelength (float) – Regular wave length, calculated from linear dispersion if set to None

  • waveType (string) – Defines regular wave theory (“Linear”, “Fenton”) Fenton: uses BCoeffs/YCoeffs provided by user

  • autoFenton (bool) –

    autoFenton=True: uses waveheight, period, depth, and g to

    calculate coeffs

    autoFenton=False: uses BCoeffs/YCoeffs provided by user

  • autoFentonOpts (dict) –

    options for autoFenton. The dictionary must contain the following entries (here the default values if autoFentonOpts is None): autoFentonOpts = {‘mode’: ‘Period’,

    ’current_criterion’: 1, ‘height_steps’: 1, ‘niter’: 40, ‘conv_crit’: 1e-05, ‘points_freesurface’: 50, ‘points_velocity’: 16, ‘points_vertical’: 20}

  • Ycoeff (numpy.ndarray) – Fenton Fourier coefficients for free-surface elevation

  • Bcoeff (numpy.ndarray) – Fenton Fourier coefficients for velocity (set to None for linear wave theory)

  • Nf (integer) – Fenton Fourier components for reconstruction (set to 1000, needs to be equal to the size of Bcoeff and Ycoeff)

  • meanVelocity (numpy.ndarray) – Mean velocity for Fenton Fourier approximation

  • phi0 (float) – Regular wave phase (0 by default)

  • fast (bool) – Switch for optimised functions

eta(x, t)[source]

Calculates free surface elevation (MonochromaticWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

mwl[source]
u(x, t)[source]

Calculates wave velocity vector (MonochromaticWaves class). :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

wavelength[source]
class proteus.WaveTools.NewWave[source]

Bases: object

This class is used for generating the NewWave theory (see Tromans et al. 1991)

Parameters
  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (numpy.ndarray) – Wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution Example for JONSWAP = {“gamma”: 3.3, “TMA”:True,”depth”: depth} TMA=True activates the TMA modification, which in turn needs the depth as a parameter

  • crestFocus (bool) – Switch to determine if crest focused or trough focused. By default true

  • xfocus (numpy array) – Position of focused crest / trough

  • tfocus (numpy array) – Time of focused crest / trough

  • Nmax (int) – Normalisation factor to get the 1/N wave event at the NewWave series

  • fast (bool) – Switch for optimised functions

Hs[source]
N[source]
Si_Jm[source]
Tlag[source]
Tp[source]
ai[source]
bandFactor[source]
depth[source]
df[source]
eta(x, t)[source]

Calculates free surface elevation (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

fi[source]
fim[source]
focus[source]
fp[source]
g[source]
gAbs[source]
kDir[source]
ki[source]
mwl[source]
omega[source]
phi[source]
tanhF[source]
tfocus[source]
u(x, t)[source]

Calculates wave velocity vector (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

vDir[source]
waveDir[source]
writeEtaSeries(Tstart, Tend, x0, fname, Lgen)[source]

Writes a timeseries of the free-surface elevation

It also returns the free surface elevation as a time-eta array. If Lgen !=[0.,0.,0.,] then Tstart is modified to account for the wave transformation at the most remote point of the relaxation zone.

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • x0 (numpy.ndarray) – Position vector of the time series

  • fname (string) – Filename for timeseries file

  • Lgen (Optional[numpy.ndarray]) – Length vector of relaxation zone

Returns

2D numpy array Nx2 containing free-surface elevation in time.

Return type

numpy.ndarray

class proteus.WaveTools.RandomWaves[source]

Bases: object

This class is used for generating plane random waves using linear reconstruction of components from a wave spectrum

Parameters
  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (numpy.ndarray) – Wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution Example for JONSWAP = {“gamma”: 3.3, “TMA”:True,”depth”: depth} TMA=True activates the TMA modification, which in turn needs the depth as a parameter

  • phi (numpy.ndarray) – Component phases (if set to None, phases are picked at random)

  • fast (bool) – Switch for optimised functions

Hs[source]
N[source]
Si_Jm[source]
Tlag[source]
Tp[source]
ai[source]
bandFactor[source]
depth[source]
df[source]
eta(x, t)[source]

Calculates free surface elevation (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

fi[source]
fim[source]
fp[source]
g[source]
gAbs[source]
kDir[source]
ki[source]
mwl[source]
omega[source]
phi[source]
tanhF[source]
u(x, t)[source]

Calculates wave velocity vector (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

vDir[source]
waveDir[source]
wavelength[source]
writeEtaSeries(Tstart, Tend, x0, fname, Lgen)[source]

Writes a timeseries of the free-surface elevation

It also returns the free surface elevation as a time-eta array. If Lgen !=[0.,0.,0.,] then Tstart is modified to account for the wave transformation at the most remote point of the relaxation zone.

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • x0 (numpy.ndarray) – Position vector of the time series

  • fname (string) – Filename for timeseries file

  • Lgen (Optional[numpy.ndarray]) – Length vector of relaxation zone

Returns

2D numpy array Nx2 containing free-surface elevation in time.

Return type

numpy.ndarray

class proteus.WaveTools.MultiSpectraRandomWaves[source]

Bases: object

This class is used for generating random waves by combining multiple spectra with different distributions and directions

Parameters
  • Nspectra (int) – Total number of spectra

  • Tp (list) – List of peak wave periods

  • Hs (list) – List of significant wave heights

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (list) – List of wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (list) – List of numbers of frequency components

  • bandFactor (list) – List of spectral band factors

  • spectName (list) – List of names of spectral distribution

  • spectral_params (list) – List of names of spectral distribution (see RandomWaves class)

  • phi (list) – List of component phases

  • fast (bool) – Switch for optimised functions

depth[source]
eta(x, t)[source]

Calculates free surface elevation (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

mwl[source]
u(x, t)[source]

Calculates wave velocity vector (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

class proteus.WaveTools.DirectionalWaves[source]

Bases: object

This class is used for generating directional random waves using linear reconstruction of components from a wave spectrum

Parameters
  • M (int) – Number of directional components

  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir0 (numpy.ndarray) – Leading wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spreadName (string) – Name of spreading distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution (see RandomWaves class)

  • spread_params (dict) – Dictionary of arguments specific to the spreading distribution Example for Cos-2s = {“s”: 10} Example for Mitsuyashu-type = {“fp”: 1/Tp, “smax”:10}

  • phi (numpy.ndarray) – Component phases (if set to None, phases are picked at random)

  • phiSymm (bool) – Switch for enabling a symmetric phase allocation across directional components

  • fast (bool) – Switch for enabling optimised functions

depth[source]
eta(x, t)[source]

Calculates free surface elevation (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

mwl[source]
u(x, t)[source]

Calculates wave velocity vector (RandomWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

class proteus.WaveTools.TimeSeries[source]

Bases: object

This class is used for generating waves from an arbirtrary free-surface elevation time series

Parameters
  • timeSeriesFile (string) – Time series file name (csv or txt)

  • skiprows (int) – Number of header rows in time series file

  • timeSeriesPosition (numpy.ndarrat) – Coordinates of the gauge / signal location

  • depth (float) – Water depth

  • N (int) – Number of frequency components

  • mwl (float) – Still water level

  • waveDir (numpy.ndarray) – Leading wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • cutoffTotal (float) – Cut off fraction, applied both at the leading and tailing parts of the series

  • rec_direct (bool) – Switch for activating direct decomposition

  • window_params (dict) – Dictionary of parameters for window method e.g. window_params = {“Nwaves”:15, “Tm”: Tp/1.1, “Window”:”costap”} (minimum parameters required) e.g. window_params = {“Nwaves”:15, “Tm”: Tp/1.1, “Window”:”costap”, “Overlap”:0.5, “Cutoff”:0.2} (full range of parameters)

  • arrayData (bool) – Switch for passing the time series as an array (False by default)

  • seriesArray (numpy.ndarray) – Free surface elevation time series given in an array format (None by default)

  • fast (bool) – Switch for enabling optimised functions

eta[source]
etaDirect(x, t)[source]

Calculates free surface elevation(Timeseries class-direct method :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

etaWindow(x, t)[source]

Calculates free surface elevation(Timeseries class-window method :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

findWindow(t)[source]

Returns the current spectral window in TimeSeries class.”

Parameters

t (float) – Time variable

Returns

Index of window as an integer

Return type

int

mwl[source]
u[source]
uDirect(x, t)[source]

Calculates wave velocity vector (Timeseries class-direct method) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

uWindow(x, t)[source]

Calculates wave velocity vector (Timeseries class-window method) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

wavelength[source]
windOut()[source]
class proteus.WaveTools.RandomWavesFast(Tstart, Tend, x0, Tp, Hs, mwl, depth, waveDir, g, N, bandFactor, spectName, spectral_params=None, phi=None, Lgen=None, Nwaves=15, Nfreq=32, checkAcc=True, fast=True)[source]

Bases: object

This class is used for generating plane random waves in an optimised manner using linear reconstruction of components from a wave spectrum

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • x0 (numpy.ndarray) – Position vector for the time series

  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (numpy.ndarray) – Wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution Example for JONSWAP = {“gamma”: 3.3, “TMA”:True,”depth”: depth} TMA=True activates the TMA modification, which in turn needs the depth as a parameter

  • phi (numpy.ndarray) – Component phases (if set to None, phases are picked at random)

  • Lgen (numpy.ndarray) – Length of the generation zone (np.array([0., 0., 0.]) by default

  • Nwaves (int) – Number of waves per window

  • Nfreq (int) – Number of Fourier components per window

  • checkAcc (bool) – Switch for enabling accuracy checks

  • fast (bool) – Switch for enabling optimised functions

printOut()[source]

Prints some properties of the time series - ONLY FOR TESTING

class proteus.WaveTools.RandomNLWaves[source]

Bases: object

This class is contains functions for calculating random waves with 2nd order corrections

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (numpy.ndarray) – Wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution Example for JONSWAP = {“gamma”: 3.3, “TMA”:True,”depth”: depth} TMA=True activates the TMA modification, which in turn needs the depth as a parameter

  • phi (numpy.ndarray) – Component phases (if set to None, phases are picked at random)

  • fast (bool) – Switch for enabling optimised functions

eta[source]
eta_2ndOrder(x, t)[source]

Calculates the free surface elevation for 2nd-order terms

Uses 2nd order random wave theory

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Returns

Free-surface elevation as a float

Return type

float

eta_linear[source]
eta_long(x, t)[source]

Calculates the free surface elevation for lower-order terms

Uses 2nd order random wave theory

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Returns

Free-surface elevation as a float

Return type

float

eta_overall(x, t, setUp)[source]

Calculates the free surface elevation with 2nd order corrections

Uses 2nd order random wave theory

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Returns

Free-surface elevation as a float

Return type

float

eta_setUp(x, t)[source]

Calculates the free surface elevation set up

Uses 2nd order random wave theory

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Returns

Free-surface elevation as a float

Return type

float

eta_short(x, t)[source]

Calculates the free surface elevation for higher-order terms

Uses 2nd order random wave theory

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Returns

Free-surface elevation as a float

Return type

float

u[source]
writeEtaSeries(Tstart, Tend, dt, x0, fname, mode, setUp, Lgen)[source]

Writes a timeseries of the free-surface elevation

It also returns the free surface elevation as a time-eta array. If Lgen !=[0.,0.,0.,] then Tstart is modified to account for the wave transformation at the most remote point of the relaxation zone.

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • dt (float) – Sampling interval

  • x0 (numpy.ndarray) – Position vector of the time series

  • fname (string) – Filename for timeseries file

  • mode (Optional[string]) – Mode of set up calculations (all, long, short, setup)

  • setUp (Optional[bool]) – Switch for activating setup calculation

  • Lgen (Optional[numpy.ndarray]) – Length vector of relaxation zone

Returns

2D numpy array Nx2 containing free-surface elevation in time.

Return type

numpy.ndarray

wtError(x, t)[source]

Raises error for using RandomNLWavesFast class instead

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

Return type

None

Raises

SystemExit

class proteus.WaveTools.RandomNLWavesFast(Tstart, Tend, x0, Tp, Hs, mwl, depth, waveDir, g, N, bandFactor, spectName, spectral_params=None, phi=None, Lgen=array([0., 0., 0.]), Nwaves=15, Nfreq=32, NLongW=10.0, fast=True)[source]

Bases: object

This class is used for generating plane random waves with 2ns order correction in an optimised manner using linear reconstruction of components from a wave spectrum

Parameters
  • Tstart (float) – Start time

  • Tend (float) – End time

  • x0 (numpy.ndarray) – Position vector for the time series

  • Tp (float) – Peak wave period

  • Hs (float) – Significant wave height

  • mwl (float) – Still water level

  • depth (float) – Water depth

  • waveDir (np.ndarray) – Wave direction vector

  • g (Numpy array) – Gravitational acceleration vector

  • N (int) – Number of frequency components

  • bandFactor (float) – Spectral band factor. fmax = bandFactor/Tp, fmin = 1/(bandFactor*Tp)

  • spectName (string) – Name of spectral distribution

  • spectral_params (dict) – Dictionary of arguments specific to the spectral distribution Example for JONSWAP = {“gamma”: 3.3, “TMA”:True,”depth”: depth} TMA=True activates the TMA modification, which in turn needs the depth as a parameter

  • phi (numpy.ndarray) – Component phases (if set to None, phases are picked at random)

  • Lgen (numpy.ndarray) – Length of the generation zone (np.array([0., 0., 0.]) by default

  • Nwaves (int) – Number of waves per window

  • Nfreq (int) – Number of Fourier components per window

  • NLongw (int) – Estmated ratio of long wave period to Tp

  • fast (bool) – Switch for enabling optimised functions

eta(x, t)[source]

Calculates free surface elevation (RandomNLWavesFast class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

u(x, t)[source]

Calculates wave velocity vector (RandomNLWavesFast class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity vector as 1D array

Return type

numpy.ndarray

class proteus.WaveTools.CombineWaves(waveList)[source]

Bases: object

This class is used for combining multiple waveTools classes, thus allowing for the generation of complex wave conditions

Parameters

waveList (list) – List of wave classes

eta(x, t)[source]

Calculates free surface elevation (combineWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Free-surface elevation as a float

Return type

float

u(x, t)[source]

Calculates wave particle velocity (combineWaves class) :param x: Position vector :type x: numpy.ndarray :param t: Time variable :type t: float

Returns

Velocity as 1D numpy array

Return type

numpy array

proteus.WaveTools.fastcos_test(phase, sinus=False)[source]

Fast cosine function with Taylor approximation - TO BE USED FOR TESTING” :param phase: Phase :type phase: double :param sinus: Switch for cosine or sine calculation :type sinus: bool

Return type

cos(phi) or sin(phi)

proteus.WaveTools.fastcosh_test(k, Z, fast=True)[source]

Fast hyperbolic cosine function with Taylor approximation - TO BE USED FOR TESTING” :param k: Wavenumber :type k: double :param Z: Z coordinate :type Z: double

Return type

cosh(k*z)

proteus.WaveTools.fastsinh_test(k, Z, fast=True)[source]

Fast hyperbolic sine function with Taylor approximation - TO BE USED FOR TESTING” :param k: Wavenumber :type k: double :param Z: Z coordinate :type Z: double

Return type

sinh(k*z)

proteus.WaveTools.coshkzd_test(k, Z, d, fast=True)[source]

Calculation of u horizontal profile cosh(k(d+Z))/sinh(kd) using fast appoximaitons and hyp trig relation cosh(a+b) = cosha*coshb+sinha*sinhb :param ———-: :param k: Wavenumber :type k: double :param Z: Z coordinate :type Z: double :param d: depth :type d: double

Return type

cosh(k*(z+d))/sinh(kd) for Z>-d/2, 0 otherwise

proteus.WaveTools.sinhkzd_test(k, Z, d, fast=True)[source]

Calculation of v vertical profile cosh(k(d+Z))/sinh(kd) using fast appoximaitons and hyp trig relation sinh(a+b) = sinha*coshb+cosha*sinhb :param ———-: :param k: Wavenumber :type k: double :param Z: Z coordinate :type Z: double :param d: depth :type d: double

Return type

sinh(k*(z+d))/sinh(kd) for Z>-d/2, 0 otherwise

proteus.WaveTools.loadExistingFunction(funcName, validFunctions)[source]

Checks if a function name is known function and returns it

Checks if a function name is present in a list of functions. If True, the function is returned. If False, raises SystemExit.

Parameters
  • funcName (string) – Function name

  • validFunctions (List[function]) – List of valid functions (list of objects)

Return type

function

Raises

SystemExit

proteus.WaveTools.setVertDir(g)[source]

Returns the unit vector for the vertical direction

The vertical direction is opposite to the gravity direction

Parameters

g (numpy.ndarray) – Gravitational acceleration vector (must have 3 components)

Return type

numpy.ndarray

proteus.WaveTools.setDirVector(vector)[source]

Returns the direction of a vector in the form of a unit vector

Parameters

vector (numpy.ndarray) – 1D numpy array with three components

Return type

numpy.ndarray

proteus.WaveTools.dirCheck(v1, v2)[source]

Checks if two vectors are vertical raises SystemError if True

Parameters
  • v1 (numpy.ndarray) – 1st vector with three components

  • v2 (numpy.ndarray) – 2nd vector with three components

Return type

None

Raises

SystemExit

proteus.WaveTools.reduceToIntervals(fi, df)[source]

Prepares the x-axis array with N elements for numerical integration

Integration is performed along he x- axis. If fi = [a1, a2, a3,…,a_N-1 a_N] then it returns the array [a1, 0.5(a1+a2), 0.5(a2+a3),…0.5(a_N-1+a_N), a_N]. Input array must have constant step

Parameters
  • fi (numpy.ndarray) – x- array with N elements

  • df (float) – Constant step of array

Return type

numpy.ndarray

proteus.WaveTools.returnRectangles(a, x)[source]

Returns 2D discrete integral array using the rectangle method

The calculation for each array element is \((\Delta y_i = 0.5(a_{n-1}+a_{n})*(x_{n-1}-x_{n})\)

Parameters
  • a (numpy.ndarray) – Description: Array of y(x) function with N+1 elements

  • x (numpy.ndarray) – Description: x- coordinate array with N elements

Return type

numpy.ndarray

proteus.WaveTools.returnRectangles3D(a, x, y)[source]

Returns 3D discrete integrals using the rectangle method

The calculation for each array element is :math: (Delta y = 0.25*(a_(n-1,m-1)+a_(n,m-1)+a_(n-1,m) … …+a_(n,m))*(x_n-1-x_n) *(z_m-1-z_m))

Parameters
  • a (numpy.ndarray) – 2D Array of y(x,y) function with (N+1)x(M+1)elements

  • x (numpy.ndarray) – Description: x- coordinate array with N+1 elements

  • y (numpy.ndarray) – Description: x- coordinate array with N+1 elements

Return type

numpy.ndarray

proteus.WaveTools.normIntegral(f, dom)[source]

Returns a normalised 2D function

The calculation is :math: (int_Omega f dOmega =1)

Parameters
  • f (numpy.ndarray) – Discrete 2D function Numpy array or list

  • dom (float) – Discrete function step

Return type

numpy.ndarray

proteus.WaveTools.eta_mode(x, t, kDir, omega, phi, amplitude)[source]

Calculates the free surface elevation for a single frequency mode

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

  • kDir (numpy.ndarray) – Wave number vector

  • omega (float) – Angular frequency

  • phi (float) – Description: Wave phase

  • amp (float) – Description: Wave amplitude

Returns

The free surface elevation at x,t

Return type

float

proteus.WaveTools.Udrift(amp, gAbs, c, d)[source]

Calculates the 2nd order Stokes drift for a linear mode

Parameters
  • amp (float) – Description: Wave amplitude

  • gAbs (float) – Magnitude of gravitational acceleration

  • c (float) – Wave celerity

  • d (float) – Water depth

Returns

Magnitude of the mean velocity drift

Return type

float

proteus.WaveTools.vel_mode(x, t, kDir, kAbs, omega, phi, amplitude, mwl, depth, vDir, gAbs)[source]

Calculates the wave velocity components for a single frequency mode

Parameters
  • x (numpy.ndarray) – Position vector

  • t (float) – Time variable

  • kDir (numpy.ndarray) – Wave number vector

  • kAbs (floatkAbs) – Wave number magnitude

  • omega (float) – Angular frequency

  • phi (float) – Description: Wave phase

  • amplidute (float) – Description: Wave amplitude

  • mwl (float) – Mean water level

  • depth (float) – Water depth

  • vDir (numpy.ndarray) – Unit vector aligned with vertical direction

Returns

1D Numpy array of the velocity vector at x,t

Return type

numpy.ndarray

proteus.WaveTools.sigma(omega, omega0)[source]

Calculates sigma function for JONSWAP spectrum

Parameters
  • omega (numpy.ndarray) – Angular frequency array

  • omega0 (numpy.ndarray) – Peak angular frequency

Returns

1D Numpy array of simga function with respect to f

Return type

numpy.ndarray

proteus.WaveTools.JONSWAP(f, f0, Hs, gamma=3.3, TMA=False, depth=None)[source]

Calculates the JONSWAP frequency spectrum (Goda 2009)

The calculation includes the TMA modification, if TMA =True

Parameters
  • f (numpy.ndarray) – Frequency array

  • f0 (float) – Peak frequency

  • Hs (float) – Significant wave height

  • gamma (Optional[float]) – Peak enhancement factor

  • TMA (bool) – Description: TMA switch

  • depth (Optional[float]) – Water depth

Returns

1D Numpy array of the spectrum in frequency domain

Return type

numpy.ndarray

proteus.WaveTools.PM_mod(f, f0, Hs)[source]

Calculates the Pierson-Moskovitz spectrum (or Bretschneider or ISSC)

Reference: http://www.orcina.com/SoftwareProducts/OrcaFlex/Documentation/Help/Content/html/Waves,WaveSpectra.htm And then to Tucker M J, 1991. Waves in Ocean Engineering. Ellis Horwood Ltd. (Chichester).

Parameters
  • f (numpy.ndarray) – Frequency array

  • f0 (float) – Peak frequency

  • Hs (float) – Significant wave height

Returns

1D Numpy array of the spectrum in frequency domain

Return type

numpy.ndarray

proteus.WaveTools.cos2s(theta, f, s=10)[source]

Calculates the cos-2s directional spreading function see USACE - CETN-I-28 http://chl.erdc.usace.army.mil/library/publications/chetn/pdf/cetn-i-28.pdf

Parameters
  • theta (numpy.ndarray) – Wave angle array

  • f (numpy.ndarray) – Frequency array

  • s (Optional[float]) – Spreading parameter

Returns

2D Numpy array of cos-2s spectrum

Return type

numpy.ndarray

proteus.WaveTools.mitsuyasu(theta, fi, f0, smax=10)[source]

The cos2s wave directional spread with wave frequency dependency

Equation from “Random Seas and Design of Maritime Structures” - Y. Goda - 2010 (3rd ed) eq. 2.22 - 2.25

Parameters
  • theta (numpy.ndarray) – Wave angle array

  • fi (numpy.ndarray) – Frequency array

  • f0 (float) – Peak frequency

  • smax (Optional[float]) – Spreading parameter

Returns

2D Numpy array of Mitsuyashu-type spectrum

Return type

numpy.ndarray

proteus.WaveTools.dispersion(w, d, g=9.81, niter=1000)[source]

Calculates the wave number for single or multiple frequencies using linear dispersion relation.

Parameters
  • w (float or np.ndarray) – Angular frequency

  • d (float) – Water depth

  • g (Optional[float]) – Gravitational acceleration

  • niter (Optional[int]) – Solution iterations

Returns

Wavenumber as a float or 1D array for multiple frequencies

Return type

float or numpy.ndarray

proteus.WaveTools.tophat(l, cutoff)[source]

Calculates and returns a top hat filter array

Parameters
  • l (int) – Length of array

  • cutoff (float) – Cut off fraction at both the leading and tailing part of the array

Return type

numpy.ndarray

proteus.WaveTools.costap(l, cutoff=0.1)[source]

Calculates and returns a top hat filter array

Parameters
  • l (int) – Length of array

  • cutoff (float) – Cut off fraction at both the leading and tailing part of the array

Return type

numpy.ndarray

proteus.WaveTools.decompose_tseries(time, eta, dt)[source]

Performs spectral analysis and calculates angular frequency components, amplitude, phase and mean level power of a time series with constant sampling.

Parameters
  • time (numpy.ndarray) – Time array

  • eta (numpy.ndarray) – Signal array

  • dt (float) – Sampling interval

Returns

A list with results with four entries:

0 -> numpy array with angular frequency components 1 -> numpy array with amplitude of each component aa 2 -> numpy array with phase of each component pp 3 -> float of the 0th fourier mode (wave setup)

Return type

List