API
Simulation Class
- class objects.simulationClass
The
simulationClass
creates asimu
object saved to the MATLAB workspace. ThesimulationClass
includes properties and methods used to define WEC-Sim’s simulation parameters and settings.- simulationClass
This method initializes the
simulationClass
.
- Constructor Summary
- Property Summary
- b2b
WEC-Sim input
- cicDt
(float) Time step to calculate Convolution Integral. Default =
dt
- cicEndTime
(float) Convolution integral time. Default =
60
s
- domainSize
(float) Size of free surface and seabed. This variable is only used for visualization. Default =
200
m
- dt
(float) Simulation time step. Default =
0.1
s
- dtOut
(float) Output sampling time. Default =
dt
- endTime
(float) Simulation end time. Default =
[]
- explorer
(string) SimMechanics Explorer ‘on’ or ‘off’. Default =
'on'
- gravity
(float) Acceleration due to gravity. Default =
9.81
m/s
- mcrMatFile
(string) mat file that contain a list of the multiple conditions runs with given conditions. Default =
[]
- mcrExcelFile
(string) File name from which to load wave statistics data. Default =
[]
- mode
(string) Simulation execution mode, ‘normal’, ‘accelerator’, ‘rapid-accelerator’. Default =
'normal'
- morisonDt
(float) Sample time to calculate Morison Element forces. Default =
dt
- nonlinearDt
(float) Sample time to calculate nonlinear forces. Default =
dt
- paraview
0 (off) , 1 (on). Default =
0
.startTime
(float) Start time for the vtk file of Paraview. Default =0
.endTime
(float) End time for the vtk file of Paraview. Default =100
.dt
(float) Timestep for Paraview. Default =0.1
.path
(string) Path of the folder for Paraview vtk files. Default ='vtk'
.- Type:
(structure) Defines the Paraview visualization.
option
(integer) Flag for paraview visualization, and writing vtp files, Options
- pressure
0 (off), 1 (on). Default =
0
- Type:
(integer) Flag to save pressure distribution, Options
- rampTime
(float) Ramp time for wave forcing. Default =
100
s
- rateTransition
‘on’, ‘off’. Default =
'on'
- Type:
(string) Flag for automatically handling rate transition for data transfer, Opyions
- reloadH5Data
0 (off), 1 (on). Default =
0
- Type:
(integer) Flag to re-load hydro data from h5 file between runs, Options
- rho
(float) Density of water. Default =
1000
kg/m^3
- saveStructure
0 (off), 1 (on). Default =
0
- Type:
(integer) Flag to save results as a MATLAB structure, Options
- saveText
0 (off), 1 (on). Default =
0
- Type:
(integer) Flag to save results as ASCII files, Options
- saveWorkspace
0 (off), 1 (on). Default =
1
- Type:
(integer) Flag to save .mat file for each run, Options
- simMechanicsFile
(string) Simulink/SimMechanics model file. Default =
'NOT DEFINED'
- solver
(string) PDE solver used by the Simulink/SimMechanics simulation. Any continuous solver in Simulink possible. Recommended to use ‘ode4, ‘ode45’ for WEC-Sim. Default =
'ode4'
- stateSpace
0 (convolution integral), 1 (state-space). Default =
0
- Type:
(integer) Flag for convolution integral or state-space calculation, Options
- FIR
0 (convolution integral), 1 (FIR). Default =
0
- Type:
(integer) Flag for FIR calculation, Options
- startTime
(float) Simulation start time. Default =
0
s
- zeroCross
(string) Disable zero cross control. Default =
'DisableAll'
- numMoorDyn
(integer) Number of moordyn blocks systems in the wec model. Default =
0
- Method Summary
Wave Class
- class objects.waveClass
The
waveClass
creates awaves
object saved to the MATLAB workspace. ThewaveClass
includes properties and methods used to define WEC-Sim’s wave input.- Constructor Summary
- Property Summary
- bem
input file
- current
0
for depth-independent model,1
for 1/7 power law variation with depth,2
for linear variation with depth, or3
for no current. Default =3
,depth
(float) Current depth [m]. Define the depth over which the sub-surface current is modeled. Must be defined for options1
and2
. The current is not calculated for any depths greater than the specified current depth. Default =0
,direction
(float) Current direction [deg]. Surface current direction defined using WEC-Sim global coordinate system. Default =0
,speed
(float) Current seed [m/s]. Surface current speed that is uniform along the water column. Default =0
- Type:
(structure) Defines the current implementation.
option
(integer) Define the sub-surface current model to be used in WEC-Sim, options include
- direction
(float) Incident wave direction(s) [deg]. Incident wave direction defined using WEC-Sim global coordinate system. Should be defined as a row vector for more than one wave direction. Default =
0
- elevationFile
(string) Data file that contains the times-series data file. Default =
'NOT DEFINED'
- gamma
(float) Defines gamma, only used for
JS
wave spectrum type. Default =[]
- height
(float) Wave height [m]. Defined as wave height for
regular
, or significant wave height forirregular
. Default ='NOT DEFINED'
- marker
Sphere,
2
: Cube,3
: Frame. Default =1
: Sphere- Type:
(structure) Defines the wave marker. loc (nx2 vector) Marker [X,Y] locations [m]. Default =
[]
.size
(float) Marker size in Pixels. Default =10
.style
Marker style, options include- Type:
1
- period
(float) Wave period [s] . Defined as wave period for
regular
, peak period forirregular
, or period of BEM data used for hydrodynamic coefficients fornoWave
. Default ='NOT DEFINED'
- phaseSeed
(integer) Defines the random phase seed, only used for
irregular
andspectrumImport
waves. Default =0
- spectrumFile
(string) Data file that contains the spectrum data file. Default =
'NOT DEFINED'
- spectrumType
PM
orJS
. Default ='NOT DEFINED'
- Type:
(string) Specifies the wave spectrum type, options include
- viz
(structure) Defines visualization options, structure contains the fields
numPointsX
for the number of visualization points in x direction, andnumPointsY
for the number of visualization points in y direction.
- waterDepth
(float) Water depth [m]. Default to BEM water depth if not set.
- spread
(float) Wave Spread probability associated with wave direction(s). Should be defined as a row vector for more than one wave direction. Default =
1
- Method Summary
- plotElevation(rampTime)
This method plots wave elevation time-history.
- Parameters:
rampTime (
float, optional
) – Specify wave ramp time to include in plot- Returns:
figure – Plot of wave elevation versus time
- Return type:
fig
- plotSpectrum()
This method plots the wave spectrum.
- Returns:
figure – Plot of wave spectrum versus wave frequency
- Return type:
fig
Body Class
- class objects.bodyClass
The
bodyClass
creates abody
object saved to the MATLAB workspace. ThebodyClass
includes properties and methods used to define WEC-Sim’s hydrodynamic and non-hydrodynamic bodies.- bodyClass
This method initializes the
bodyClass
and creates abody
object.- Parameters:
h5File (
string
) – String specifying the location of the body h5 file- Returns:
body – bodyClass object
- Return type:
obj
- Constructor Summary
- Property Summary
- adjMassFactor
WEC-Sim input
- centerBuoyancy
(3x1 float vector) Body center of buoyancy [m]. Defined in the following format [x y z]. For hydrodynamic bodies this is defined in the h5 file while for nonhydrodynamic bodies this is defined by the user. Default =
[]
.
- centerGravity
(3x1 float vector) Body center of gravity [m]. Defined in the following format [x y z]. For hydrodynamic bodies this is defined in the h5 file while for nonhydrodynamic bodies this is defined by the user. Default =
[]
.
- dof
(integer) Number of degree of freedoms (DOFs). For hydrodynamic bodies this is given in the h5 file. If not defined in the h5 file, Default =
6
.
- excitationIRF
(vector) Defines excitation Impulse Response Function, only used with the waveClass
elevationImport
type. Default =[]
.
- flex
0 (off) or 1 (on). Default =
0
.- Type:
(integer) Flag for flexible body, Options
- gbmDOF
(integer) Number of degree of freedoms (DOFs) for generalized body mode (GBM). Default =
[]
.
- geometryFile
(string) Path to the body geometry
.stl
file.
- h5File
(char array, string, cell array of char arrays, or cell array or strings) hdf5 file containing the hydrodynamic data
- hydroStiffness
(6x6 float matrix) Linear hydrostatic stiffness matrix. If the variable is nonzero, the matrix will override the h5 file values. Create a 3D matrix (6x6xn) for variable hydrodynamics. Default =
zeros(6)
.
- inertia
(1x3 float vector) Rotational inertia or mass moment of inertia [kg*m^{2}]. Defined by the user in the following format [Ixx Iyy Izz]. Default =
[]
.
- inertiaProducts
(1x3 float vector) Rotational inertia or mass products of inertia [kg*m^{2}]. Defined by the user in the following format [Ixy Ixz Iyz]. Default =
[]
.
- initial
(structure) Defines the initial displacement of the body.
displacement
(1x3 float vector) is defined as the initial displacement of the body center of gravity (COG) [m] in the following format [x y z], Default = [0 0 0
].axis
(1x3 float vector) is defined as the axis of rotation in the following format [x y z], Default = [0 1 0
].angle
(float) is defined as the initial angular displacement of the body COG [rad], Default =0
.
- linearDamping
(6x6 float matrix) Defines linear damping coefficient matrix. Create a 3D matrix (6x6xn) for variable hydrodynamics. Default =
zeros(6)
.
- mass
(float) Translational inertia or mass [kg]. Defined by the user or specify ‘equilibrium’ to set the mass equal to the fluid density times displaced volume. Default =
[]
.
- meanDrift
0 (no), 1 (yes, from control surface) or 2 (yes, from momentum conservation). Default =
0
.- Type:
(integer) Flag for mean drift force, Options
- morisonElement
0 (off), 1 (on) or 2 (on), Default =
0
, Option 1 uses an approach that allows the user to define drag and inertial coefficients along the x-, y-, and z-axes and Option 2 uses an approach that defines the Morison Element with normal and tangential tangential drag and interial coefficients.cd
(1x3 float vector) is defined as the viscous normal and tangential drag coefficients in the following format, Option 1[cd_x cd_y cd_z]
, Option 2[cd_N cd_T NaN]
, Default =[0 0 0]
.ca
is defined as the added mass coefficent for the Morison Element in the following format, Option 1[ca_x ca_y ca_z]
, Option 2[ca_N ca_T NaN]
, Default =[0 0 0]
,area
is defined as the characteristic area for the Morison Element [m^2] in the following format, Option 1[Area_x Area_y Area_z]
, Option 2[Area_N Area_T NaN]
, Default =[0 0 0]
.VME
is the characteristic volume of the Morison Element [m^3], Default =0
.rgME
is defined as the vector from the body COG to point of application for the Morison Element [m] in the following format[x y z]
, Default =[0 0 0]
.z
is defined as the unit normal vector center axis of the Morison Element in the following format, Option 1 not used, Option 2[n_{x} n_{y} n_{z}]
, Default =[0 0 1]
.- Type:
(structure) Defines the Morison Element properties connected to the body.
option
(integer) for Morison Element calculation, Options
- name
(string) Specifies the body name. For hydrodynamic bodies this is defined in h5 file. For nonhydrodynamic bodies this is defined by the user, Default =
[]
.
- nonHydro
0 (hydro body), 1 (non-hydro body), 2 (drag body). Default =
0
.- Type:
(integer) Flag for non-hydro body, Options
- nonlinearHydro
0 (linear), 1 (nonlinear), 2 (nonlinear). Default =
0
- Type:
(integer) Flag for nonlinear hydrohanamics calculation, Options
- quadDrag
(structure) Defines the viscous quadratic drag forces. Create an array of structures for variable hydrodynamics. First option define
drag
, (6x6 float matrix), Default =zeros(6)
. Second option definecd
, (1x6 float vector), Default =[0 0 0 0 0 0]
, andarea
, (1x6 float vector), Default =[0 0 0 0 0 0]
.
- QTFs
0 (off), 1 (on). Default =
0
- Type:
(integer) Flag for QTFs calculation, Options
- paraview
0 (no) or 1 (yes). Default =
1
, only called in paraview.- Type:
(integer) Flag for visualisation in Paraview either, Options
- variableHydro
(structure) Defines the variable hydro implementation.
option
(float) Flag to turn variable hydrodynamics on or off.hydroForceIndexInitial
(float) Defines the initial value of the hydroForceIndex, which should correspond to the hydroForce data (cg, cb, volume, water depth, valid cicEndTime, added mass integrated with the body during runtime) and h5File of the body at equilibrium.
- viz
(structure) Defines visualization properties in either SimScape or Paraview.
color
(1x3 float vector) is defined as the body visualization color, Default = [1 1 0
].opacity
(integer) is defined as the body opacity, Default =1
.
- volume
(float) Displaced volume at equilibrium position [m^{3}]. For hydrodynamic bodies this is defined in the h5 file while for nonhydrodynamic bodies this is defined by the user. Default =
[]
.
- yaw
0 (off), 1 (on). Default =
0
.threshold
(float) Yaw position threshold (in degrees) above which excitation coefficients will be interpolated in passive yaw. Default =1
[deg].- Type:
(structure) Defines the passive yaw implementation.
option
(integer) Flag for passive yaw calculation, Options
- Method Summary
- setInitDisp(relCoord, axisAngleList, addLinDisp)
Function to set a body’s initial displacement
This function assumes that all rotations are about the same relative coordinate. If not, the user should input a relative coordinate of 0,0,0 and use the additional linear displacement parameter to set the centerGravity or location correctly.
- Parameters:
relCoord (
[1 3] float vector
) – Distance from x_rot to the body center of gravity or the constraint or pto location as defined by: relCoord = centerGravity - x_rot. [m]axisAngleList (
[nAngle 4] float vector
) – List of axes and angles of the rotations with the format: [n_x n_y n_z angle] (angle in rad) Rotations applied consecutively in order of dimension 1addLinDisp (
[1 3] float vector
) – Initial linear displacement (in addition to the displacement caused by rotation) [m]
- plotStl()
Method to plot the body .stl mesh and normal vectors.
- calculateForceAddedMass(acc)
This method calculates and outputs the real added mass force time history. This encompasses both the contributions of the added mass coefficients and applied during simulation, and the component from added mass that is lumped with the body mass during simulation.
This function must be called after body.restoreMassMatrix() and body.storeForceAddedMass()
- Parameters:
obj (
bodyClass
) – Body whose added mass force is being updatedacc (
float array
) – Timeseries of the acceleration at each simulation time step
- Returns:
actualAddedMassForce – Time series of the actual added mass force
- Return type:
float array
Constraint Class
- class objects.constraintClass
The
constraintClass
creates aconstraint
object saved to the MATLAB workspace. TheconstraintClass
includes properties and methods used to define constraints between the body motion relative to the global reference frame or relative to other bodies.- constraintClass
This method initilizes the
constraintClass
and creates aconstraint
object.- Parameters:
filename (
string
) – String specifying the name of the constraint- Returns:
constraint – contraintClass object
- Return type:
obj
- Constructor Summary
- Property Summary
- hardStops
input file
- initial
(structure) Defines the initial displacement of the constraint.
displacement
(3x1 float vector) is defined as the initial displacement of the constraint [m] in the following format [x y z], Default = [0 0 0
].
- extension
(string) Specify priority level for extension. `` ‘High’ or ‘Low’ `` Default =
High
.
- location
(1x3 float vector) Constraint location [m]. Defined in the following format [x y z]. Default =
[999 999 999]
.
- name
(string) Specifies the constraint name. For constraints this is defined by the user, Default =
NOT DEFINED
.
- orientation
(structure) Defines the orientation axis of the constraint.
z
(1x3 float vector) defines the direciton of the Z-coordinate of the constraint, Default = [0 0 1
].y
(1x3 float vector) defines the direciton of the Y-coordinate of the constraint, Default = [0 1 0
].x
(3x1 float vector) internally calculated vector defining the direction of the X-coordinate for the constraint, Default =[]
.rotationMatrix
(3x3 float matrix) internally calculated rotation matrix to go from standard coordinate orientation to the constraint coordinate orientation, Default =[]
.
- Method Summary
- checkInputs()
This method checks WEC-Sim user inputs and generates error messages if parameters are not properly defined.
- setInitDisp(relCoord, axisAngleList, addLinDisp)
Function to set a constraints’s initial displacement
This function assumes that all rotations are about the same relative coordinate. If not, the user should input a relative coordinate of 0,0,0 and use the additional linear displacement parameter to set the cg or location correctly.
- Parameters:
relCoord (
[1 3] float vector
) – Distance from x_rot to the body center of gravity or the constraint or pto location as defined by: relCoord = cg - x_rot. [m]axisAngleList (
[nAngle 4] float vector
) – List of axes and angles of the rotations with the format: [n_x n_y n_z angle] (angle in rad) Rotations applied consecutively in order of dimension 1addLinDisp (
[1 3] float vector
) – Initial linear displacement (in addition to the displacement caused by rotation) [m]
PTO Class
- class objects.ptoClass
The
ptoClass
creates apto
object saved to the MATLAB workspace. TheptoClass
includes properties and methods used to define PTO connections between the body motion relative to the global reference frame or relative to other bodies.- ptoClass
This method initilizes the
ptoClass
and creates apto
object.- Parameters:
filename (
string
) – String specifying the name of the pto- Returns:
pto – ptoClass object
- Return type:
obj
- Constructor Summary
- Property Summary
- damping
input file
- equilibriumPosition
(float) Linear PTO equilibrium position, m or deg. Default = 0.
- hardStops
(float) Define lower limit transition region, over which spring and damping values ramp up to full values. Increase for stability. m or deg.
Default = 1e-4
- initial
(structure) Defines the initial displacement of the pto.
displacement
(3x1 float vector) is defined as the initial displacement of the pto [m] in the following format [x y z], Default = [0 0 0
].
- extension
(string) Specify priority level for extension. `` ‘High’ or ‘Low’ `` Default =
High
.
- location
(3x1 float vector) PTO location [m]. Defined in the following format [x y z]. Default =
[999 999 999]
.
- name
(string) Specifies the pto name. For ptos this is defined by the user, Default =
NOT DEFINED
.
- orientation
(structure) Defines the orientation axis of the pto.
z
(1x3 float vector) defines the direction of the Z-coordinate of the pto, Default = [0 0 1
].y
(1x3 float vector) defines the direction of the Y-coordinate of the pto, Default = [0 1 0
].x
(1x3 float vector) internally calculated vector defining the direction of the X-coordinate for the pto, Default =[]
.rotationMatrix
(3x3 float matrix) internally calculated rotation matrix to go from standard coordinate orientation to the pto coordinate orientation, Default =[]
.
- pretension
(float) Linear PTO pretension, N or N-m. Default = 0.
- stiffness
(float) Linear PTO stiffness coefficient. Default = 0.
- Method Summary
- checkInputs()
This method checks WEC-Sim user inputs and generates error messages if parameters are not properly defined.
- setInitDisp(relCoord, axisAngleList, addLinDisp)
Function to set a pto’s initial displacement
This function assumes that all rotations are about the same relative coordinate. If not, the user should input a relative coordinate of 0,0,0 and use the additional linear displacement parameter to set the cg or location correctly.
- Parameters:
relCoord (
[1 3] float vector
) – Distance from x_rot to the body center of gravity or the constraint or pto location as defined by: relCoord = cg - x_rot. [m]axisAngleList (
[nAngle 4] float vector
) – List of axes and angles of the rotations with the format: [n_x n_y n_z angle] (angle in rad) Rotations applied consecutively in order of dimension 1addLinDisp (
[1 3] float vector
) – Initial linear displacement (in addition to the displacement caused by rotation) [m]
Mooring Class
- class objects.mooringClass
The
mooringClass
creates amooring
object saved to the MATLAB workspace. ThemooringClass
includes properties and methods used to define cable connections relative to other bodies. It is suggested that themooringClass
be used for connections between bodies and the global reference frame.This class contains mooring parameters and settings
- Constructor Summary
- Property Summary
- initial
input file
- location
(1x3 float vector) Mooring Reference location. Default =
[0 0 0]
- matrix
(structure) Defines the mooring parameters.
damping
(6x6 float matrix) Matrix of damping coefficients, Default =zeros(6)
.stiffness
(6x6 float matrix) Matrix of stiffness coefficients, Default =zeros(6)
.preTension
(1x6 float vector) Array of pretension force in each dof, Default =[0 0 0 0 0 0]
.
- moorDyn
(integer) Flag to indicate and initialize a MoorDyn connection, 0 or 1. Default =
0
- moorDynLines
(integer) Number of lines in MoorDyn. Default =
0
- moorDynNodes
(integer) number of nodes for each line. Default =
'NOT DEFINED'
- name
(string) Name of the mooring. Default =
'NOT DEFINED'
- moorDynInputFile
(string) Name of the MoorDyn input file path. Outputs will be written to this path. Default =
Mooring/lines.txt
- lookupTableFlag
(integer) Flag to indicate a mooring look-up table, 0 or 1. Default =
0
- lookupTableFile
(string) Mooring look-up table file name. Default =
''
;
- lookupTable
(array) Lookup table data. Force data in 6 DOFs indexed by displacement in 6 DOF
- Method Summary
- checkInputs()
This method checks WEC-Sim user inputs and generates error messages if parameters are not properly defined.
- setInitDisp(relCoord, axisAngleList, addLinDisp)
Function to set a mooring’s initial displacement
This function assumes that all rotations are about the same relative coordinate. If not, the user should input a relative coordinate of 0,0,0 and use the additional linear displacement parameter to set the cg or orientation correctly.
- Parameters:
relCoord (
[1 3] float vector
) – Distance from x_rot to the body center of gravity or the constraint or pto location as defined by: relCoord = cg - x_rot. [m]axisAngleList (
[nAngle 4] float vector
) – List of axes and angles of the rotations with the format: [n_x n_y n_z angle] (angle in rad) Rotations applied consecutively in order of dimension 1addLinDisp (
[1 3] float vector
) – Initial linear displacement (in addition to the displacement caused by rotation) [m]
- loadLookupTable()
Method to load the lookup table and assign to the mooringClass
- callMoorDynLib()
Initialize MoorDyn Lib (Windows:dll or OSX:dylib)
- closeMoorDynLib()
Close MoorDyn Lib
Cable Class
- class objects.cableClass
The
cableClass
creates acable
object saved to the MATLAB workspace. ThecableClass
includes properties and methods used to define cable connections relative to other bodies. It is suggested that thecableClass
be used for connections between joints or ptos.- cableClass
This method initilizes the
cableClass
and creates acable
object.- Parameters:
name (
string
) – String specifying the name of the cablebaseConnection (
string
) – Variable for the base constraint/pto as a stringfollowerConnection (
string
) – Variable for the follower constraint/pto as a string
- Returns:
cable – cableClass object
- Return type:
obj
- Constructor Summary
- Property Summary
- damping
input file
- inertia
(1x3 float vector) body moments of inertia kg-m^2, default [1 1 1]
- inertiaProducts
(1x3 float vector) body products of inertia kg-m^2, default [0 0 0]
- initial
(structure) Defines the initial displacement of the body.
displacement
(3x1 float vector) is defined as the initial displacement of the body center of gravity (COG) [m] in the following format [x y z], Default = [0 0 0
].axis
(3x1 float vector) is defined as the axis of rotation in the following format [x y z], Default = [0 1 0
].angle
(float) is defined as the initial angular displacement of the body COG [rad], Default =0
.
- cableLength
(float) Cable equilibrium length (m), calculated from rotloc and preTension. Default =`0`.
- linearDamping
(1x6 float vector) linear damping aplied to body motions
- mass
(float) mass in kg, default 1
- name
(string) Defines the Cable name. Default =
NOT DEFINED
.
- orientation
(structure) Defines the orientation axis of the pto.
z
(1x3 float vector) defines the direciton of the Z-coordinate of the pto, Default = [0 0 1
].y
(1x3 float vector) defines the direciton of the Y-coordinate of the pto, Default = [0 1 0
].x
(1x3 float vector) internally calculated vector defining the direction of the X-coordinate for the pto, Default =[]
.rotationMatrix
(3x3 float matrix) internally calculated rotation matrix to go from standard coordinate orientation to the pto coordinate orientation, Default =[]
.
- paraview
(integer) Flag for visualisation in Paraview either 0 (no) or 1 (yes). Default =
1
since only called in paraview.
- preTension
(float) Cable pretension (N).
- quadDrag
(structure) Defines the viscous quadratic drag forces. First option define
drag
, (6x6 float matrix), Default =zeros(6)
. Second option definecd
, (6x1 float vector), Default =zeros(6,1)
, andarea
, (6x1 float vector), Default =zeros(6,1)
.
- stiffness
(float) Cable stiffness (N/m). Default = 0.
- viz
(structure) Defines visualization properties in either SimScape or Paraview.
color
(3x1 float vector) is defined as the body visualization color, Default = [1 1 0
].opacity
(integer) is defined as the body opacity, Default =1
.
- base
internal
- follower
(structure) Defines the follower connection. centerBuoyancy (1 x 3 float vector) center of buoyancy location of the base drag body, Default = [0 0 0]. centerGravity (1 x 3 float vector) center of gravity location of the base drag body, Default = [0 0 0]. location (3x1 float vector) base location [m], Defined in the following format [x y z], Default =
[]
. name (string) name of the base constraint or PTO, Default = ‘NOT DEFINED’;
- Method Summary
- checkInputs()
This method checks WEC-Sim user inputs and generates error messages if parameters are not properly defined.
- setInitDisp(relCoord, axisAngleList, addLinDisp)
Function to set a body’s initial displacement
This function assumes that all rotations are about the same relative coordinate. If not, the user should input a relative coordinate of 0,0,0 and use the additional linear displacement parameter to set the cg or location correctly.
- Parameters:
relCoord (
[1 3] float vector
) – Distance from x_rot to the body center of gravity or the constraint or pto location as defined by: relCoord = cg - x_rot. [m]axisAngleList (
[nAngle 4] float vector
) – List of axes and angles of the rotations with the format: [n_x n_y n_z angle] (angle in rad) Rotations applied consecutively in order of dimension 1addLinDisp (
[1 3] float vector
) – Initial linear displacement (in addition to the displacement caused by rotation) [m]
- setCg()
This method specifies the Cg of the drag bodies as coincident with the fixed ends of the cable, if not otherwise specied.
Response Class
- class objects.responseClass
The
responseClass
creates anoutput
object saved to the MATLAB workspace that contains structures for each instance of a WEC-Sim class (e.g.waveClass
,bodyClass
,constraintClass
,ptoClass
,cableClass
,mooringClass
, etc).- responseClass
This method initializes the
responseClass
, reads output from each instance of a WEC-Sim class (e.g.waveClass
,bodyClass
,ptoClass
,mooringClass
, etc) and saves the response to anoutput
object.- Returns:
output – responseClass object
- Return type:
obj
- wave
This property contains a structure for each instance of the
waveClass
type
(string) = ‘waveType’time
(array) = [# of time-steps x 1]elevation
(array) = [# of time-steps x 1]waveGauge1Elevation
(array) = [# of time-steps x 1] Wave elevation at the location of wave gauge 1waveGauge2Elevation
(array) = [# of time-steps x 1] Wave elevation at the location of wave gauge 2waveGauge3Elevation
(array) = [# of time-steps x 1] Wave elevation at the location of wave gauge 3
- bodies
This property contains a structure for each instance of the
bodyClass
(i.e. for each Body block)name
(string) = ‘bodyName’time
(array) = [# of time-steps x 1]position
(array) = [# of time-steps x 6]velocity
(array) = [# of time-steps x 6]acceleration
(array) = [# of time-steps x 6]forceTotal
(array) = [# of time-steps x 6] The sum of all hydrodynamic forces applied to the bodyforceExcitation
(array) = [# of time-steps x 6] The sum of the Froude-Krylov excitation force and the mean drift force exerted by the incoming wave on the bodyforceRadiationDamping
(array) = [# of time-steps x 6] The negative radiation damping force due to body velocityforceAddedMass
(array) = [# of time-steps x 6] The negative added mass force due to body accelerationforceRestoring
(array) = [# of time-steps x 6] The negative sum of the gravity force, buoyant force, the hydrostatic stiffness force, and any moment due to separation between the center of gravity and the center of buoyancyforceMorisonAndViscous
(array) = [# of time-steps x 6] The negative sum of the Morison element force and the viscous drag forceforceLinearDamping
(array) = [# of time-steps x 6] The negative force due to linear damping and the body velocity
There are 4 additional
output.bodies
arrays when using nonlinear hydro and Paraview output:cellPressures_time
(array) = [# of Paraview time-steps x 1] Nonlinear calculation timeseriescellPressures_hydrostatic
(array) = [# of Paraview time-steps x # of mesh faces] Hydrostatic pressure on each stl facetcellPressures_waveLinear
(array) = [# of Paraview time-steps x # of mesh faces] Excitation pressure on each stl facet given zero displacement and the mean free surfacecellPressures_waveNonLinear
(array) = [# of Paraview time-steps x # of mesh faces] Excitation pressure on each stl facet given the instantaneous displacement and instantaneous free surface
- constraints
This property contains a structure for each instance of the
constraintClass
(i.e. for each Constraint block). Constraint motion is relative from frame F to frame B. Constraint forces act on frame F.name
(string) = ‘constraintName’time
(array) = [# of time-steps x 1]position
(array) = [# of time-steps x 6] The constraint position relative to the initial conditionvelocity
(array) = [# of time-steps x 6] The constraint velocity relative to the initial conditionacceleration
(array) = [# of time-steps x 6] The constraint acceleration relative to the initial conditionforceConstraint
(array) = [# of time-steps x 6] The force required to resist motion in the restricted DOFs
- ptos
This property contains a structure for each instance of the
ptoClass
(i.e. for each PTO block). PTO motion is relative from frame F to frame B. PTO forces act on frame F.name
(string) = ‘ptoName’time
(array) = [# of time-steps x 1]position
(array) = [# of time-steps x 6] The constraint position relative to the initial conditionvelocity
(array) = [# of time-steps x 6] The constraint velocity relative to the initial conditionacceleration
(array) = [# of time-steps x 6] The constraint acceleration relative to the initial conditionforceTotal
(array) = [# of time-steps x 6] The sum of the actuation, constraint and internal mechanics forcesforceActuation
(array) = [# of time-steps x 6] The prescribed force input to the PTO to control its motionforceConstraint
(array) = [# of time-steps x 6] The force required to resist motion in the restricted DOFsforceInternalMechanics
(array) = [# of time-steps x 6] The net force in the joint DOF due to stiffness and dampingpowerInternalMechanics
(array) = [# of time-steps x 6] The net power lost in the joint DOF due to stiffness and damping
- cables
This property contains a structure for each instance of the
cableClass
(i.e. for each Cable block)name
(string) = ‘cableName’time
(array) = [# of time-steps x 1]position
(array) = [# of time-steps x 6]velocity
(array) = [# of time-steps x 6]acceleration
(array) = [# of time-steps x 6]forceTotal
(array) = [# of time-steps x 6] The sum of the actuation and constraint forcesforceActuation
(array) = [# of time-steps x 6] The cable tensionforceConstraint
(array) = [# of time-steps x 6] The force required to resist motion in the restricted DOFs
- mooring
This property contains a structure for each instance of the
mooringClass
using the mooring matrix (i.e. for each MooringMatrix block)position
(array) = [# of time-steps x 6]velocity
(array) = [# of time-steps x 6]forceMooring
(array) = [# of time-steps x 6] The sum of the stiffness, damping and pretension forces applied on the body by the mooring
- moorDyn
This property contains a structure for each instance of the
mooringClass
using MoorDyn (i.e. for each MoorDyn block)Lines
(struct) = [1 x 1] Contains the time and fairlead tensionsLine#
(struct) = [1 x 1] One structure for each mooring line: Line1, Line2, etc. Each line structure contains node positions in x, y, z and segment tensions
- ptoSim
This property contains a structure for each instance of the
ptoSimClass
(i.e. for each PTO-Sim block).time
(struct) = [# of time-steps x 1] Simulation timeseries
There are additional
output.ptoSim
structs corresponding to the Simulink blocks used:pistonCF
(struct) = [1 x # of pistons] Structure containing timeseries of compressible fluid piston properties including absolute power, force, position, velocitypistonNCF
(array) = [1 x # of pistons] Structure containing timeseries of non-compressible fluid piston properties including absolute power, force, top pressure and bottom pressurecheckValve
(struct) = [1 x # of valves] Structure containing timeseries of check valve properies including volume flow ratevalve
(struct) = [1 x # of valves] Structure containing timeseries of valve properties including volume flow rateaccumulator
(struct) = [1 x # of accumulators] Structure containing timeseries of accumulator properties including pressure and volumehydraulicMotor
(struct) = [1 x # of motors] Structure containing timeseries of hydraulic motor properties including angular velocity and volume flow raterotaryGenerator
(struct) = [1 x # of generators] Structure containing timeseries of rotary generator properties including electrical power and generated powersimpleDD
(struct) = [1 x # of generators] Structure containing timeseries of direct drive PTO properties including forces and electrical powerpmLinearGenerator
(struct) = [1 x # of generators] Structure containing timeseries of permanent magnet linear generator properties including absolute power, force, friction force, current, voltage, velocity and electrical powerpmRotaryGenerator
(struct) = [1 x # of generators] Structure containing timeseries of permanent magnet rotary generator properties including absolute power, force, friction force, current, voltage, velocity and electrical powermotionMechanism
(struct) = [1 x # of mechanisms] Structure containing timeseries of motion mechanism properties including PTO torque, angular position and angular velocity
- windTurbine
This property contains a structure for each instance of the
windTurbineClass
name
(string) = ‘windTurbineName’time
(array) = [# of time-steps x 1]windSpeed
(array) = [# of time-steps x 1]turbinePower
(array) = [# of time-steps x 1]rotorSpeed
(array) = [# of time-steps x 1]bladePitch
(array) = [# of time-steps x 1] Pitch position of blade 1nacelleAcceleration
(array) = [# of time-steps x 1]towerBaseLoad
(array) = [# of time-steps x 6] 6DOF force at the constraint between the floating body and tower basetowerTopLoad
(array) = [# of time-steps x 6] 6DOF force at the constraint between the tower base and tower nacellebladeRootLoad
(array) = [# of time-steps x 1] 6DOF force at the constraint between blade 1 and the hubgenTorque
(array) = [# of time-steps x 1] Torque on the generatorazimuth
(array) = [# of time-steps x 1] azimuthal angle of the generator
- Constructor Summary
- Property Summary
- ptoSim
This property contains a structure for each instance of the
ptoSimClass
(i.e. for each PTO-Sim block).
- windTurbine
This property contains a structure for each instance of the
windTurbineClass
- Method Summary
- plotResponse(bodyNum, comp)
This method plots the response of a body for a given DOF.
- Parameters:
bodyNum (
integer
) – the body number to plotcomp (
integer
) – the response component (i.e. dof) to be plotted (e.g. 1-6)
- plotForces(bodyNum, comp)
This method plots the forces on a body for a given DOF.
- Parameters:
bodyNum (
integer
) – the body number to plotcomp (
integer
) – the force component (i.e. dof) to be plotted (e.g. 1-6)
- saveViz(simu, body, waves, options)
This method plots and saves the wave elevation and body geometry over time to visualize the waves and response
- Parameters:
simu – simulationClass object
body – bodyClass object
waves – waveClass object
options –
- axisLimits1x6 float matrix
x, y, and z-bounds of figure axes in the format: [min x, max x, min y, max y, min z, max z] Default = [-simu.domainSize/2 simu.domainSize/2 -simu.domainSize/2 simu.domainSize/2 -waves.waterDepth maximumHeight]
- timesPerFramefloat
number of simulation time steps per video frame (higher number decreases computation time) Default = 1
- startEndTime1x2 float matrix
Array defining the start and end times of the visualization Default = [min(t) max(t)]
- saveSettinginteger
0 = video, 1 = gif. Default = 0
- writeText()
This method writes WEC-Sim outputs to a (ASCII) text file. This method is executed by specifying
simu.outputtxt=1
in thewecSimInputFile.m
.