API Documentation
Provide core functionality for solving the pseudo-spectral problem.
- class wecopttool.core.WEC(fb: FloatingBody, mass: ndarray, hydrostatic_stiffness: ndarray, f0: float, nfreq: int, dissipation: Optional[ndarray] = None, stiffness: Optional[ndarray] = None, f_add: Optional[Mapping[str, Callable[[WEC, ndarray, ndarray], ndarray]]] = None, constraints: list[dict] = [], rho: float = 1025.0, depth: float = inf, g: float = 9.81)[source]
Class representing a specific wave energy converter (WEC). An instance contains the following information about the WEC, environment, and dynamic modeling:
Geometry
Degrees of freedom
Mass properties
Hydrostatic properties
Linear frequency domain hydrodynamic coefficients
Water properties
Additional dynamic forces (power take-off, mooring, nonlinear hydrodynamics, etc.)
Constraints
- __init__(fb: FloatingBody, mass: ndarray, hydrostatic_stiffness: ndarray, f0: float, nfreq: int, dissipation: Optional[ndarray] = None, stiffness: Optional[ndarray] = None, f_add: Optional[Mapping[str, Callable[[WEC, ndarray, ndarray], ndarray]]] = None, constraints: list[dict] = [], rho: float = 1025.0, depth: float = inf, g: float = 9.81) None [source]
- Parameters
fb (capytaine.FloatingBody) – The WEC as a capytaine floating body (mesh + DOFs).
mass (np.ndarray) – Mass matrix shape of (
ndof
xndof
).hydrostatic_stiffness (np.ndarray) – Hydrostatic stiffness matrix matrix of shape (
ndof
xndof
).f0 (float) – Initial frequency (in Hz) for frequency array. Frequency array given as [
f0
, 2f0
, …,nfreq f0
].nfreq (int) – Number of frequencies in frequency array. See
f0
.dissipation (np.ndarray) – Additional dissipation for the impedance calculation in
capytaine.post_pro.impedance
. Shape: (ndof
xndof
x1
) or (ndof
xndof
xnfreq
).stiffness (np.ndarray) – Additional stiffness for the impedance calculation in
capytaine.post_pro.impedance
. Shape: (ndof
xndof
x1
) or (ndof
xndof
xnfreq
).f_add (dict[str, Callable]) – Additional forcing terms (e.g. buoyancy, gravity, PTO, mooring, etc.) for the WEC dynamics in the time-domain. Dictionary entries should be
entry = {'name': function_handle}
. Takes three inputs: (1) the WEC object, (2) the WEC dynamics state (1D np.ndarray), and (3) the optimization state (1D np.ndarray) and outputs the force time-series (1D np.ndarray).constraints (list[dict]) – Constraints for the constrained optimization. See
scipy.optimize.minimize
.rho (float, optional) – Water density in \(kg/m^3\). Default is 1025.
depth (float, optional) – Water depth in \(m\). Default is np.inf
g (float, optional) – Gravitational acceleration in \(m/s^2\). Default is 9.81.
- __setattr__(name, value)[source]
Delete dependent attributes when user manually modifies an attribute.
- bem_calc_impedance(tol=1e-06)[source]
Calculate the impedance, ensure positive real diagonal, and create impedance MIMO matrix.
- decompose_decision_var(state: ndarray) tuple[ndarray, ndarray] [source]
Split the state vector into the WEC dynamics state and the optimization (control) state. x = [x_wec, x_opt].
- property derivative_mat
Derivative matrix for the state vector.
- dofmat_to_vec(mat: ndarray) ndarray [source]
Flatten a matrix that has one column per DOF. Opposite of
vec_to_dofmat
.
- property dt
Time spacing.
- property f0
Initial frequency (and spacing) in Hz. See
freq
.
- property f_add
Additional forces on the WEC (e.g., PTO, mooring, buoyancy, gravity)
- property freq
Frequency array f=[f0, 2*f0, …, nfreq*f0] in Hz.
- static from_file(fpath: str | Path) WEC [source]
Create a WEC instance from a file saved using WEC.to_file.
- initial_x_wec_guess_analytic(waves: xr.Dataset) np.ndaray [source]
Initial guess for x_wec based on optimal hydrodynamic solution to be passed to wec.solve.
- Parameters
waves (xr.Dataset) – Wave DataSet
- Returns
Initial guess for x_wec
- Return type
x_wec_0
Examples
>>> x_wec_0 = wec.initial_x_wec_guess_analytic(regular_wave) >>> wec.solve(regular_wave, obj_fun=pto.average_power, nstate_opt=pto.nstate, scale_x_wec=1.0, scale_x_opt=0.01, scale_obj=1e-1, x_wec_0=x_wec_0, )
- make_time_mat(nsubsteps: int = 1, include_mean: bool = True) ndarray [source]
Assemble the time matrix that converts the state to time-series.
- Parameters
nsubsteps (int) – Number of subdivisions between the default (implied) time steps.
include_mean (bool) – Whether the state vector includes a mean component.
- Returns
time_mat
- Return type
np.ndarray
- make_time_vec(nsubsteps: int = 1) ndarray [source]
Assemble the time vector with n subdivisions.
- Parameters
nsubsteps (int) – Number of subdivisions between the default (implied) time steps.
- Returns
time_vec
- Return type
np.ndarray
- natural_frequency() tuple[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], int] [source]
Return natural frequency or frequencies.
See wot.natural_frequency().
- property ncomponents
Number of state values for each DOF in the WEC dynamics.
- property ndof
Number of degrees of freedom of the WEC.
- property nfreq
Number of frequencies in frequency array. See
freq
.
- property nstate_wec
Length of the WEC dynamics state vector.
- property omega
Frequency array in radians per second ω=2πf.
- optimal_position(waves: xr.DataSet) np.ndarray [source]
Return optimal position spectrum for hydrodynamic problem.
See wot.optimal_position()
- optimal_velocity(waves: xr.DataSet) np.ndarray [source]
Return optimal velocity spectrum for hydrodynamic problem.
See wot.optimal_velocity()
- property period
math:T=1/f in seconds.
- Type
Period
- plot_impedance(style: str = 'Bode', option: str = 'symmetric', show: bool = True)[source]
Plot impedance.
See wot.plot_impedance().
- power_limit(waves: xr.DataSet) np.ndarray [source]
Return theoretical power limit for hydrodynamic problem.
See wot.power_limit()
- read_bem(fpath: str | Path, tol: float = 1e-06) None [source]
Read a BEM solution from a NetCDF file.
- Parameters
fpath (str) – Name of file to read BEM data from.
tol (float) – Minimum value for the diagonal terms of (radiation damping + dissipation).
- run_bem(wave_dirs: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]] = [0], tol: float = 1e-06) None [source]
Run the BEM for the specified wave directions.
See
wot.run_bem
.- Parameters
wave_dirs (list[float]) – List of wave directions to evaluate BEM at (degrees).
tol (float) – Minimum value for the diagonal terms of (radiation damping + dissipation).
- solve(waves: xr.Dataset, obj_fun: Callable[[WEC, np.ndarray, np.ndarray], float], nstate_opt: int, x_wec_0: Optional[np.ndarray] = None, x_opt_0: Optional[np.ndarray] = None, scale_x_wec: Optional[list] = None, scale_x_opt: npt.ArrayLike | float = 1.0, scale_obj: float = 1.0, optim_options: dict[str, Any] = {}, use_grad: bool = True, maximize: bool = False, bounds_wec: Optional[Bounds] = None, bounds_opt: Optional[Bounds] = None, unconstrained_first: Optional[bool] = False, callback: Callable[[np.ndarray]] = None) tuple[xr.Dataset, xr.Dataset, np.ndarray, np.ndarray, float, OptimizeResult] [source]
Solve the WEC co-design problem.
- Parameters
waves (xr.Dataset) – The wave, described by two 2D DataArrays: elevation variance S (m^2*s) and phase phase (radians) with coordinates of radial frequency omega (radians) and wave direction wave_direction (radians). The frequencies and wave directions must match those in the bem_data in self.hydro.
obj_fun (function) – Objective function for the control optimization. Takes three inputs: (1) the WEC object, (2) the WEC dynamics state (1D np.ndarray), and (3) the optimization state (1D np.ndarray) and outputs the scalar objective function: tuple[WEC, np.ndarray, np.ndarray] -> float.
nstate_opt (int) – Length of the optimization (controls) state vector.
x_wec_0 (np.ndarray) – Initial guess for the WEC dynamics state. If
None
it is randomly initiated.x_opt_0 (np.ndarray) – Initial guess for the optimization (control) state. If
None
it is randomly initiated.scale_x_wec (list) – Factors to scale each DOF in
x_wec
by, to improve convergence. List lengthndof
.scale_x_opt (npt.ArrayLike | float) – Factor to scale
x_opt
by, to improve convergence. A single float or an array of sizenstate_opt
.scale_obj (float) – Factor to scale
obj_fun
by, to improve convergence.optim_options (dict) – Optimization options passed to the optimizer. See
scipy.optimize.minimize
.use_grad (bool) – Whether to use gradient information in the optimization.
maximize (bool) – Whether to maximize the objective function. The default is
False
to minimize the objective function.bounds_wec (Bounds) – Bounds on the WEC components of the decsision variable; see scipy.optimize.minimize
bounds_opt (Bounds) – Bounds on the optimization (control) components of the decsision variable; see scipy.optimize.minimize
unconstrained_first (bool) – If True, run
solve
without constraints to get scaling and initial guess. The default is False.callback (function) – Called after each iteration; see scipy.optimize.minimize. The default is reported via logging at the INFO level.
- Returns
time_dom (xr.Dataset) – Dataset containing the time-domain results.
freq_dom (xr.Dataset) – Dataset containing the frequency-domain results.
x_wec (np.ndarray) – Optimal WEC state.
x_opt (np.ndarray) – Optimal control state.
objective (float) – optimized value of the objective function.
res (optimize.optimize.OptimizeResult) – Raw optimization results.
- td_to_fd(td: ndarray) ndarray [source]
Convert from frequency domain to time domain using the iFFT.
- property tf
Final time (period).
- property time
Time array.
- property time_mat
Matrix to convert from the state vector to a time-series.
- vec_to_dofmat(vec: ndarray) ndarray [source]
Convert a vector back to a matrix with one column per DOF. Opposite of
dofmat_to_vec
.
- property w0
Initial frequency (and spacing) in rad/s. See
freq
.
- wecopttool.core.complex_to_real_amplitudes(fd: ndarray) ndarray [source]
Convert from one complex amplitude to two real amplitudes per frequency.
- wecopttool.core.complex_xarray_from_netcdf(fpath: str | Path) xr.Dataset [source]
Read a NetCDF file with complex entries as an xarray DataSet.
- wecopttool.core.complex_xarray_to_netcdf(fpath: str | Path, bem_data: xr.Dataset) None [source]
Save an xarray dataSet with complex entries as a NetCDF file.
- wecopttool.core.freq_array(f0: float, nfreq: int) ndarray [source]
Construct equally spaced frequency array.
- wecopttool.core.natural_frequency(impedance: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], freq: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]]) tuple[Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], int] [source]
Find the natural frequency based on the lowest magnitude impedance.
- Parameters
impedance (np.ndarray) – Complex impedance matrix. Shape:
nfreq
xndof
xndof
freq (list[float]) – Frequencies.
- Returns
f_n (float) – Natural frequency.
ind (int) – Index of natural frequency.
- wecopttool.core.optimal_position(excitation: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], impedance: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], omega: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]]) ndarray [source]
Find optimal position.
- Parameters
excitation (np.ndarray) – Complex excitation spectrum. Shape:
nfreq
xndof
impedance (np.ndarray) – Complex impedance matrix. Shape:
nfreq
xndof
xndof
- Returns
Optimal position for power absorption.
- Return type
optimal_position
- wecopttool.core.optimal_velocity(excitation: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], impedance: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]]) ndarray [source]
Find optimal velocity.
- Parameters
excitation (np.ndarray) – Complex excitation spectrum. Shape:
nfreq
xndof
impedance (np.ndarray) – Complex impedance matrix. Shape:
nfreq
xndof
xndof
- Returns
Optimal velocity for power absorption.
- Return type
opt_vel
- wecopttool.core.plot_impedance(impedance: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], freq: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], style: str = 'Bode', option: str = 'diagonal', show: bool = False, dof_names: Optional[list[str]] = None) tuple[Figure, ndarray] [source]
Plot the impedance matrix.
- Parameters
impedance (np.ndarray) – Complex impedance matrix. Shape:
nfreq
xndof
xndof
freq (list[float]) – Frequencies in Hz.
style ({'Bode','complex'}) – Whether to plot magnitude and angle (
Bode
) or real and imaginary (complex
) parts.option ({'diagonal', 'symmetric', 'all'}) – Which terms of the matrix to plot: ‘diagonal’ to plot only the diagonal terms, ‘symmetric’ to plot only the lower triangular terms, and ‘all’ to plot all terms.
show (bool) – Whether to show the figure.
dof_names (list[str]) –
- Returns
fig (matplotlib.figure.Figure)
axs (np.ndarray[matplotlib.axes._subplots.AxesSubplot])
- wecopttool.core.post_process_continuous_time(results: DataArray) Callable[[float], float] [source]
Create a continuous function from the results in an xarray DataArray.
The DataArray must be indexed by “omega”: frequency in rad/s. There should be no other indices.
- Parameters
results (xr.DataArray) – DataArray containing the pseudo-spectral results.
- Returns
func – Continuous-time function.
- Return type
Callable
- wecopttool.core.power_limit(excitation: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], impedance: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]]) ndarray [source]
Find upper limit for power.
- Parameters
excitation (np.ndarray) – Complex excitation spectrum. Shape:
nfreq
xndof
impedance (np.ndarray) – Complex impedance matrix. Shape:
nfreq
xndof
xndof
- Returns
Upper limit for power absorption.
- Return type
power_limit
- wecopttool.core.real_to_complex_amplitudes(fd: ndarray, first_row_is_mean: bool = True) ndarray [source]
Convert from two real amplitudes to one complex amplitude per frequency.
- wecopttool.core.run_bem(fb: FloatingBody, freq: Iterable[float] = [inf], wave_dirs: Iterable[float] = [0], rho: float = 1025.0, g: float = 9.81, depth: float = inf, write_info: Iterable[str] = []) Dataset [source]
Run Capytaine for a range of frequencies and wave directions.
- Parameters
fb (capytaine.FloatingBody) – The WEC as a Capytaine floating body (mesh + DOFs).
freq (list[float]) – List of frequencies to evaluate BEM at.
wave_dirs (list[float]) – List of wave directions to evaluate BEM at (degrees).
rho (float, optional) – Water density in \(kg/m^3\).
g (float, optional) – Gravitational acceleration in \(m/s^2\).
depth (float, optional) – Water depth in \(m\).
write_info (list[str], optional) – List of information to keep, passed to capytaine solver. Options are: wavenumber, wavelength, mesh, hydrostatics.
- Returns
BEM results from capytaine.
- Return type
xarray.Dataset
- wecopttool.core.scale_dofs(scale_list: list[float], ncomponents: int) ndarray [source]
Create a scaling vector based on a different scale for each DOF.
- Parameters
scale_list (list) – Scale for each DOF.
ncomponents (int) – Number of elements in the state vector for each DOF.
- Returns
np.ndarray
- Return type
Scaling vector.
- wecopttool.core.wave_excitation(bem_data: Dataset, waves: Dataset, time_mat: ndarray) tuple[Dataset, Dataset] [source]
Compute the frequency- and time-domain wave excitation force.
- Parameters
bem_data (xarray.Dataset) – BEM data for the WEC obtained from capytaine.
waves (xarray.Dataset) – The wave, described by two 2D DataArrays: elevation variance S (m^2*s) and phase phase (radians) with coordinates of radial frequency omega (radians) and wave direction wave_direction (radians). The frequencies and wave directions must match those in the bem_data.
- Returns
freq_dom (xarray.Dataset) – Frequency domain wave excitation and elevation.
time_dom (xarray.Dataset) – Time domain wave excitation and elevation.