MoorDyn
Public Member Functions | Public Attributes | List of all members
moordyn::Line Class Referencefinal

A mooring line. More...

#include <Line.hpp>

Inheritance diagram for moordyn::Line:
Inheritance graph
[legend]
Collaboration diagram for moordyn::Line:
Collaboration graph
[legend]

Public Member Functions

 Line (moordyn::Log *log, size_t lineId)
 Constructor. More...
 
 ~Line ()
 Destructor.
 
void setup (int number, LineProps *props, real l, unsigned int n, EnvCondRef env_in, shared_ptr< ofstream > outfile, string channels, real dtM0)
 Setup a line. More...
 
void setWaves (moordyn::WavesRef waves_in, moordyn::SeafloorRef seafloor_in)
 Set the environmental data. More...
 
void initialize ()
 Initialize the line object. More...
 
void initialize (InstanceStateVarView state)
 Sets the initial line state. More...
 
void setInitialTmax (moordyn::real Tmax0)
 Set the initial Tmax values. More...
 
void setInitialTmean (moordyn::real Tmean0)
 Set the initial mean tension values. More...
 
elastic_model getElasticModel () const
 Get the elastic model used by the line. More...
 
unsigned int getN () const
 Number of segments. More...
 
moordyn::real getUnstretchedLength () const
 Get the unstretched length of the line. More...
 
void setUnstretchedLength (const moordyn::real len)
 Set the unstretched length of the line. More...
 
void setUnstretchedLengthVel (const moordyn::real v)
 Set the unstretched length rate of change of the line. More...
 
void updateUnstretchedLength (const moordyn::real dt=0.0)
 Update the unstretched length of the line, according to the velocity. More...
 
bool isConstantEA () const
 Get whether the line is governed by a non-linear stiffness or a constant one. More...
 
moordyn::real getConstantEA () const
 Get the constant stiffness of the line. More...
 
void setConstantEA (moordyn::real EA_in)
 Set the constant stiffness of the line. More...
 
const vecgetNodePos (unsigned int i) const
 Get the position of a node. More...
 
const vecgetNodeVel (unsigned int i) const
 Get the velocity of a node. More...
 
const vecgetNodeForce (unsigned int i) const
 Get the net force on a node. More...
 
const vec getNodeTen (unsigned int i) const
 Get the tension on a node, including the internal line damping. More...
 
const vecgetNodeBendStiff (unsigned int i) const
 Get the tension on a node, including the internal line damping. More...
 
const vecgetNodeWeight (unsigned int i) const
 Get the weight and bouyancy force acting on the node. More...
 
const vec getNodeDrag (unsigned int i) const
 Get the drag force acting on the node. More...
 
const vec getNodeFroudeKrilov (unsigned int i) const
 Get the Froude-Krilov force acting on the node. More...
 
const vec getNodeSeabedForce (unsigned int i) const
 Get the sea bed reaction force acting on the node. More...
 
const realgetNodeCurv (unsigned int i) const
 Get the line curvature at a node position. More...
 
const matgetNodeM (unsigned int i) const
 Get the mass and added mass matrix. More...
 
std::vector< vecgetNodeCoordinates () const
 Get the array of coordinates of all nodes along the line. More...
 
void getFASTtens (float *FairHTen, float *FairVTen, float *AnchHTen, float *AnchVTen) const
 Get the tensions at the fairlead and anchor in a FASTv7 friendly way. More...
 
void getEndStuff (vec &Fnet_out, vec &Moment_out, mat &M_out, EndPoints end_point) const
 Get the force, moment and mass at the line endpoint. More...
 
real GetLineOutput (OutChanProps outChan)
 Get line output. More...
 
void storeWaterKin (real dt, std::vector< std::vector< moordyn::real >> zeta, std::vector< std::vector< moordyn::real >> f, std::vector< std::vector< vec >> u, std::vector< std::vector< vec >> ud)
 store wave/current kinematics time series for this line More...
 
real calcSubSeg (unsigned int firstNodeIdx, unsigned int secondNodeIdx, real surfaceHeight)
 calculate the volume of the segment between firstNodeIdx and secondNodeIdx submerged More...
 
std::pair< real, realgetDrag () const
 Get the drag coefficients. More...
 
void setDrag (real cdn, real cdt)
 Set the drag coefficients. More...
 
void scaleDrag (real scaler)
 Scale the drag coefficients. More...
 
void enablePb ()
 Enable the pressure bending forces (disabled by default) More...
 
void disablePb ()
 Disable the pressure bending forces (disabled by default)
 
bool enabledPb () const
 Check if pressure bending forces are considered. More...
 
void setPin (std::vector< real > p)
 Set the internal pressure at each node of the line. More...
 
void setTime (real time)
 Set the line simulation time. More...
 
- Public Member Functions inherited from moordyn::Instance
 Instance (moordyn::Log *log)
 Costructor. More...
 
virtual ~Instance ()=default
 Destructor.
 
const size_t id () const
 Get the unique identifier of this instance. More...
 
- Public Member Functions inherited from moordyn::io::IO
 IO (moordyn::Log *log)
 Costructor. More...
 
virtual ~IO ()=default
 Destructor.
 
void Save (const std::string filepath)
 Save the entity into a file. More...
 
void Load (const std::string filepath)
 Loads the entity from a file. More...
 
- Public Member Functions inherited from moordyn::LogUser
 LogUser (Log *log=NULL)
 Constructor. More...
 
 ~LogUser ()
 Destructor.
 
void SetLogger (Log *log)
 Set the log handler. More...
 
LogGetLogger () const
 Get the log handler. More...
 

Public Attributes

int number
 Line ID.
 
size_t lineId
 
bool IC_gen = false
 boolean to indicate dynamic relaxation occuring
 
void setState (const InstanceStateVarView r)
 Set the line state. More...
 
void setEndKinematics (vec r, vec rd, EndPoints end_point)
 Set the position and velocity of an end point. More...
 
void setEndOrientation (vec q, EndPoints end_point, EndPoints rod_end_point)
 set end node unit vector More...
 
vec getEndSegmentMoment (EndPoints end_point, EndPoints rod_end_point) const
 Get the bending moment at the end point. More...
 
void getStateDeriv (InstanceStateVarView drdt)
 Calculate forces and get the derivative of the line's states. More...
 
void Output (real)
 
const size_t stateN () const
 Get the number of state variables required by this instance. More...
 
const size_t stateDims () const
 Get the dimension of the state variable. More...
 
std::vector< uint64_t > Serialize (void)
 Produce the packed data to be saved. More...
 
uint64_t * Deserialize (const uint64_t *data)
 Unpack the data to restore the Serialized information. More...
 
void saveVTK (const char *filename)
 Save the line on a VTK (.vtp) file. More...
 
const leanvtk::VTPWriter * getVTK () const
 Get the VTK writer. More...
 

Additional Inherited Members

- Protected Member Functions inherited from moordyn::io::IO
ofstream MakeFile (const std::string filepath) const
 Create an output file and write the MoorDyn magic header. More...
 
std::tuple< uint64_t, uint64_t * > LoadFile (const std::string filepath) const
 Open an input file and load the data. More...
 
uint64_t Serialize (const uint64_t &i)
 Pack an unsigned integer to make it writable. More...
 
uint64_t Serialize (const int64_t &i)
 Pack an integer to make it writable. More...
 
uint64_t Serialize (const real &f)
 Pack a float to make it writable. More...
 
std::vector< uint64_t > Serialize (const vec &m)
 Pack a 3D vector to make it writable. More...
 
std::vector< uint64_t > Serialize (const vec6 &m)
 Pack a 6D vector to make it writable. More...
 
std::vector< uint64_t > Serialize (const mat &m)
 Pack a 3x3 matrix to make it writable. More...
 
std::vector< uint64_t > Serialize (const mat6 &m)
 Pack a 6x6 matrix to make it writable. More...
 
std::vector< uint64_t > Serialize (const quaternion &m)
 Pack a quaternion to make it writable. More...
 
std::vector< uint64_t > Serialize (const XYZQuat &m)
 Pack an XYZQuat to make it writable. More...
 
std::vector< uint64_t > Serialize (const std::vector< real > &l)
 Pack a list of floating point numbers to make it writable. More...
 
std::vector< uint64_t > Serialize (const std::vector< vec > &l)
 Pack a list of 3D vectors to make it writable. More...
 
std::vector< uint64_t > Serialize (const std::vector< vec6 > &l)
 Pack a list of 6D vectors to make it writable. More...
 
std::vector< uint64_t > Serialize (const std::vector< mat > &l)
 Pack a list of 3x3 matrices to make it writable. More...
 
std::vector< uint64_t > Serialize (const std::vector< mat6 > &l)
 Pack a list of 6x6 matrices to make it writable. More...
 
std::vector< uint64_t > Serialize (const Eigen::Matrix< real, Eigen::Dynamic, Eigen::Dynamic > &l)
 Pack an arbitrarily large matrix. More...
 
template<typename T >
std::vector< uint64_t > Serialize (const std::vector< std::vector< T >> &l)
 Pack a list of lists to make it writable This function might act recursively. More...
 
uint64_t * Deserialize (const uint64_t *in, uint64_t &out)
 Unpack a loaded unsigned integer. More...
 
uint64_t * Deserialize (const uint64_t *in, int64_t &out)
 Unpack a loaded integer. More...
 
uint64_t * Deserialize (const uint64_t *in, real &out)
 Unpack a loaded floating point number. More...
 
uint64_t * Deserialize (const uint64_t *in, vec &out)
 Unpack a loaded 3D vector. More...
 
uint64_t * Deserialize (const uint64_t *in, vec6 &out)
 Unpack a loaded 6D vector. More...
 
uint64_t * Deserialize (const uint64_t *in, mat &out)
 Unpack a loaded 3x3 matrix. More...
 
uint64_t * Deserialize (const uint64_t *in, mat6 &out)
 Unpack a loaded 6x6 matrix. More...
 
uint64_t * Deserialize (const uint64_t *in, quaternion &out)
 Unpack a loaded quaternion. More...
 
uint64_t * Deserialize (const uint64_t *in, XYZQuat &out)
 Unpack a loaded XYZQuat. More...
 
uint64_t * Deserialize (const uint64_t *in, std::vector< real > &out)
 Unpack a loaded list of floating point numbers. More...
 
uint64_t * Deserialize (const uint64_t *in, std::vector< vec > &out)
 Unpack a loaded list of 3D vectors. More...
 
uint64_t * Deserialize (const uint64_t *in, std::vector< vec6 > &out)
 Unpack a loaded list of 6D vectors. More...
 
uint64_t * Deserialize (const uint64_t *in, std::vector< mat > &out)
 Unpack a loaded list of 3x3 matrices. More...
 
uint64_t * Deserialize (const uint64_t *in, std::vector< mat6 > &out)
 Unpack a loaded list of 6x6 matrices. More...
 
uint64_t * Deserialize (const uint64_t *in, Eigen::Matrix< real, Eigen::Dynamic, Eigen::Dynamic > &out)
 Unpack an arbitrarily large matrix. More...
 
template<typename T >
uint64_t * Deserialize (const uint64_t *in, std::vector< std::vector< T >> &out)
 Unpack a loaded list of lists. More...
 
- Protected Attributes inherited from moordyn::LogUser
Log_log
 The log handler.
 

Detailed Description

A mooring line.

Each mooring line is divided on a set of nodes connected by segments, with two points at the end points

[point (node 0)] - seg 0 - [node 1] - ... - seg n-1 - [point (node N)]

Depending on the length of the line, $ l $, the number of segments, $ n $, and the material density and Young's modulus, $ \rho $ & $ E $, a natural oscillation period can be defined:

$ T = \frac{l}{\pi n} \sqrt{\frac{\rho}{E}} $

The integration time step (moordyn::MoorDyn.dtM0) should be smaller than this natural period to avoid numerical instabilities

Constructor & Destructor Documentation

◆ Line()

moordyn::Line::Line ( moordyn::Log log,
size_t  lineId 
)

Constructor.

Parameters
logLogging handler
lineIdthe incremental Id of the line

Member Function Documentation

◆ calcSubSeg()

real moordyn::Line::calcSubSeg ( unsigned int  firstNodeIdx,
unsigned int  secondNodeIdx,
real  surfaceHeight 
)

calculate the volume of the segment between firstNodeIdx and secondNodeIdx submerged

This must be used with adjacent nodes for accurate results. It is currently only implemented for the still water case.

Parameters
firstNodeIdxIndex of the first node of the segment
secondNodeIdxIndex of the second node of the segment
surfaceHeightHeight of the water surface (assumed locally horizontal)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Deserialize()

uint64_t * moordyn::Line::Deserialize ( const uint64_t *  data)
virtual

Unpack the data to restore the Serialized information.

This is the inverse of Serialize(void)

Parameters
dataThe packed data
Returns
A pointer to the end of the file, for debugging purposes

Implements moordyn::io::IO.

Here is the call graph for this function:

◆ enabledPb()

bool moordyn::Line::enabledPb ( ) const
inline

Check if pressure bending forces are considered.

Returns
true if pressure bending forces are considered, false otherwise

◆ enablePb()

void moordyn::Line::enablePb ( )
inline

Enable the pressure bending forces (disabled by default)

The internal pressure can be provided at each node by calling ::internalPress(), while the external pressure is computed as the hydrostatic pressure plus the dynamic pressure (see moordyn::Waves::getWaveKinLine()).

If no internal pressure is provided, zeros will be considered.

◆ getConstantEA()

moordyn::real moordyn::Line::getConstantEA ( ) const
inline

Get the constant stiffness of the line.

This value is useless if non-linear stiffness is considered

Returns
The static stiffness EA value
See also
::IsConstantEA()

◆ getDrag()

std::pair<real, real> moordyn::Line::getDrag ( ) const
inline

Get the drag coefficients.

Returns
The normal (transversal) and tangential (axial) drag coefficients

◆ getElasticModel()

elastic_model moordyn::Line::getElasticModel ( ) const
inline

Get the elastic model used by the line.

Returns
The elastic model. See ::elastic_model

◆ getEndSegmentMoment()

vec moordyn::Line::getEndSegmentMoment ( EndPoints  end_point,
EndPoints  rod_end_point 
) const

Get the bending moment at the end point.

The bending moment is defined as $\bar{r} \frac{E I}{\vert \bar{r} \vert^2}$, with $\bar{r}$ the vector joining the endpoint and the next line node

This method is already taking into account the line and rod end points to appropriately set the sign of the resulting moment

Parameters
end_pointEither ENDPOINT_B or ENDPOINT_A
rod_end_pointEither ENDPOINT_B or ENDPOINT_A
Returns
The moment vector
Exceptions
invalid_value_errorIf either end_point or rod_end_point are not valid end point qualifiers
Here is the call graph for this function:

◆ getEndStuff()

void moordyn::Line::getEndStuff ( vec Fnet_out,
vec Moment_out,
mat M_out,
EndPoints  end_point 
) const
inline

Get the force, moment and mass at the line endpoint.

Parameters
Fnet_outThe output force vector
Moment_outThe output moment vector
M_outThe output mass matrix
end_pointEither ENDPOINT_TOP or ENDPOINT_BOTTOM
Exceptions
invalid_value_errorIf end_point is not a valid end point qualifier

◆ getFASTtens()

void moordyn::Line::getFASTtens ( float *  FairHTen,
float *  FairVTen,
float *  AnchHTen,
float *  AnchVTen 
) const
inline

Get the tensions at the fairlead and anchor in a FASTv7 friendly way.

Parameters
FairHTenHorizontal tension on the fairlead
FairVTenVertical tension on the fairlead
AnchHTenHorizontal tension on the anchor
AnchVTenVertical tension on the anchor

◆ GetLineOutput()

real moordyn::Line::GetLineOutput ( OutChanProps  outChan)

Get line output.

This function is useful when outputs are set in the line properties

Parameters
outChanThe output channel/field
Returns
The output value, 0.0 if a non-valid field is set
Here is the call graph for this function:

◆ getN()

unsigned int moordyn::Line::getN ( ) const
inline

Number of segments.

The number of nodes can be computed as moordyn::Line::getN() + 1

Returns
The number of segments, moordyn::Line::N
Here is the caller graph for this function:

◆ getNodeBendStiff()

const vec& moordyn::Line::getNodeBendStiff ( unsigned int  i) const
inline

Get the tension on a node, including the internal line damping.

If it is an inner node, the average of the tension at the surrounding segments is provided. If the node is a line-end, the associated ending segment tension is provided

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeCoordinates()

std::vector<vec> moordyn::Line::getNodeCoordinates ( ) const
inline

Get the array of coordinates of all nodes along the line.

Returns
The positions array

◆ getNodeCurv()

const real& moordyn::Line::getNodeCurv ( unsigned int  i) const
inline

Get the line curvature at a node position.

Parameters
iThe line node index
Returns
The line curvature (1 / m)
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1
Note
The curvature is only computed if bending stiffness (moordyn::Line::EI) is not zero. Otherwise the curvature of every single node will be zero.

◆ getNodeDrag()

const vec moordyn::Line::getNodeDrag ( unsigned int  i) const
inline

Get the drag force acting on the node.

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeForce()

const vec& moordyn::Line::getNodeForce ( unsigned int  i) const
inline

Get the net force on a node.

The net force is the total force acting over a line node.

To get the different components of the force use ::getNodeTen() , ::getNodeBendStiff() , ::getNodeWeight() , ::getNodeDrag() , ::getNodeFroudeKrilov() and ::getNodeSeaBedForce()

Note
The net force is not the sum of all the components that you cat extract from the API. For instance, the tension contribution on the internal nodes is the difference between the tensions of the adjacent segments, while ::getNodeTen() is returning the average.
Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1
Here is the caller graph for this function:

◆ getNodeFroudeKrilov()

const vec moordyn::Line::getNodeFroudeKrilov ( unsigned int  i) const
inline

Get the Froude-Krilov force acting on the node.

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeM()

const mat& moordyn::Line::getNodeM ( unsigned int  i) const
inline

Get the mass and added mass matrix.

Parameters
iThe line node index
Returns
The mass matrix
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodePos()

const vec& moordyn::Line::getNodePos ( unsigned int  i) const
inline

Get the position of a node.

Parameters
iThe line node index
Returns
The position
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeSeabedForce()

const vec moordyn::Line::getNodeSeabedForce ( unsigned int  i) const
inline

Get the sea bed reaction force acting on the node.

This is eventually including the friction force

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeTen()

const vec moordyn::Line::getNodeTen ( unsigned int  i) const
inline

Get the tension on a node, including the internal line damping.

If it is an inner node, the average of the tension at the surrounding segments is provided. If the node is a line-end, the associated ending segment tension is provided. This does not account for the direction of pull. In the calculation of Fnet for the node, the direction of pull is accounted for.

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1
Here is the caller graph for this function:

◆ getNodeVel()

const vec& moordyn::Line::getNodeVel ( unsigned int  i) const
inline

Get the velocity of a node.

Parameters
iThe line node index
Returns
The velocity
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getNodeWeight()

const vec& moordyn::Line::getNodeWeight ( unsigned int  i) const
inline

Get the weight and bouyancy force acting on the node.

Parameters
iThe line node index
Returns
The tension
Exceptions
invalid_value_errorIf the node index i is bigger than the number of nodes, moordyn::Line::N + 1

◆ getStateDeriv()

void moordyn::Line::getStateDeriv ( InstanceStateVarView  drdt)
virtual

Calculate forces and get the derivative of the line's states.

Parameters
drdtOutput state derivatives
Exceptions
nan_errorIf nan values are detected in any node position

Implements moordyn::Instance.

Here is the call graph for this function:

◆ getUnstretchedLength()

moordyn::real moordyn::Line::getUnstretchedLength ( ) const
inline

Get the unstretched length of the line.

Returns
The unstretched length, moordyn::Line::UnstrLen

◆ getVTK()

const leanvtk::VTPWriter* moordyn::Line::getVTK ( ) const
inline

Get the VTK writer.

This function is useful for writing multiblock .vtm files

Returns
The VTK .vtu writer

◆ initialize() [1/2]

void moordyn::Line::initialize ( )

Initialize the line object.

Note
Becasue Line.hpp is the only function that calls this, no state vectors need to be returned. This is different from the structure of the other objects initialize functions.
Exceptions
moordyn::output_file_errorIf an outfile has been provided, but it cannot be written
invalid_value_errorIf there is no enough water depth
Here is the call graph for this function:

◆ initialize() [2/2]

void moordyn::Line::initialize ( InstanceStateVarView  state)
inlinevirtual

Sets the initial line state.

Parameters
stateThe line state for initializing (see Line::setState for the state structure0)
Note
This calls Line::Initialize()

Implements moordyn::Instance.

◆ isConstantEA()

bool moordyn::Line::isConstantEA ( ) const
inline

Get whether the line is governed by a non-linear stiffness or a constant one.

Returns
true if the stiffness of the line is constant, false if a non-linear stiffness has been set

◆ saveVTK()

void moordyn::Line::saveVTK ( const char *  filename)

Save the line on a VTK (.vtp) file.

Parameters
filenameThe output file name
Exceptions
output_file_errorIf If the file cannot be saved
Here is the call graph for this function:

◆ scaleDrag()

void moordyn::Line::scaleDrag ( real  scaler)
inline

Scale the drag coefficients.

Parameters
scalerThe drag coefficients scale factor

◆ Serialize()

std::vector< uint64_t > moordyn::Line::Serialize ( void  )
virtual

Produce the packed data to be saved.

The produced data can be used afterwards to restore the saved information afterwards calling Deserialize(void).

Thus, this function is not processing the information that is extracted from the definition file

Returns
The packed data

Implements moordyn::io::IO.

Here is the call graph for this function:

◆ setConstantEA()

void moordyn::Line::setConstantEA ( moordyn::real  EA_in)
inline

Set the constant stiffness of the line.

This value is useless if non-linear stiffness is considered

Parameters
EAThe constant axial stiffness EA value (N)
See also
::IsConstantEA()

◆ setDrag()

void moordyn::Line::setDrag ( real  cdn,
real  cdt 
)
inline

Set the drag coefficients.

Parameters
cdnNormal (transversal) coefficient
cdttangential (axial) coefficient

◆ setEndKinematics()

void moordyn::Line::setEndKinematics ( vec  r,
vec  rd,
EndPoints  end_point 
)

Set the position and velocity of an end point.

Parameters
rPosition
rdVelocity
end_pointEither ENDPOINT_TOP or ENDPOINT_BOTTOM
Exceptions
invalid_value_errorIf end_point is not a valid end point qualifier

◆ setEndOrientation()

void moordyn::Line::setEndOrientation ( vec  q,
EndPoints  end_point,
EndPoints  rod_end_point 
)

set end node unit vector

This method is called by an eventually attached Rod, only applicable for bending stiffness

Parameters
qThe direction unit vector
end_pointEither ENDPOINT_B or ENDPOINT_A
rod_end_pointEither ENDPOINT_B or ENDPOINT_A
Exceptions
invalid_value_errorIf either end_point or rod_end_point are not valid end point qualifiers

◆ setInitialTmax()

void moordyn::Line::setInitialTmax ( moordyn::real  Tmax0)
inline

Set the initial Tmax values.

Parameters
Tmax0Initial Tmax value

◆ setInitialTmean()

void moordyn::Line::setInitialTmean ( moordyn::real  Tmean0)
inline

Set the initial mean tension values.

Parameters
Tmean0Initial Tmean value

◆ setPin()

void moordyn::Line::setPin ( std::vector< real p)

Set the internal pressure at each node of the line.

If this is not provided, the last stored values (zero by default) will be considered.

This function is useless unless ::enablePb() is called.

Parameters
pThe pressure values (Pa) at each node
Exceptions
invalid_value_errorIf p has wrong size

◆ setState()

void moordyn::Line::setState ( const InstanceStateVarView  r)
virtual

Set the line state.

Parameters
rThe line state matrix. See Line::setState in Line.cpp for structure
Note
End node kinematics are not handled here
See also
moordyn::Line::setEndState

Implements moordyn::Instance.

◆ setTime()

void moordyn::Line::setTime ( real  time)
inline

Set the line simulation time.

Parameters
timeSimulation time

◆ setUnstretchedLength()

void moordyn::Line::setUnstretchedLength ( const moordyn::real  len)
inline

Set the unstretched length of the line.

This function is consistently changing the unstretched length of each segment, moordyn::Line::l

Parameters
lenThe unstretched length, moordyn::Line::UnstrLen
Note
This function should be called after moordyn::Line::initialize()
Warning
The lines damping is not changed, which might affect the stability

◆ setUnstretchedLengthVel()

void moordyn::Line::setUnstretchedLengthVel ( const moordyn::real  v)
inline

Set the unstretched length rate of change of the line.

Parameters
vThe unstretched length rate of change, moordyn::Line::UnstrLend
See also
moordyn::Line::setUnstretchedLength()

◆ setup()

void moordyn::Line::setup ( int  number,
LineProps props,
real  l,
unsigned int  n,
EnvCondRef  env_in,
shared_ptr< ofstream >  outfile,
string  channels,
real  dtM0 
)

Setup a line.

Parameters
numberLine number
propsLine properties
lUnstretched line length
nNumber of segments
env_inGlobal struct that holds environmental settings
outfileThe outfile where information shall be written
channelsThe channels/fields that shall be printed in the file
dtM0moordyn internal timestep. Used in VIV
Here is the caller graph for this function:

◆ setWaves()

void moordyn::Line::setWaves ( moordyn::WavesRef  waves_in,
moordyn::SeafloorRef  seafloor_in 
)
inline

Set the environmental data.

Parameters
waves_inGlobal Waves object
seafloor_inGlobal 3D Seafloor object

◆ stateDims()

const size_t moordyn::Line::stateDims ( ) const
inlinevirtual

Get the dimension of the state variable.

Returns
3 components for positions and 3 components for velocities, i.e. 6 components. An additional component is returned if the VIV model or viscoelastic model is activated.
Note
If VIV and vsicoelastic, return 8 if VIV xor viscoelastic, return 7 if normal, return 6 See comments in Line::setState to see the line state structure
Warning
This function shall be called after ::setup()

Reimplemented from moordyn::Instance.

◆ stateN()

const size_t moordyn::Line::stateN ( ) const
inlinevirtual

Get the number of state variables required by this instance.

Returns
The number of internal nodes if no viscoelastic model otherwise, return the number of segments
Note
See comments in Line::setState to see the line state structure
Warning
This function shall be called after ::setup()

Reimplemented from moordyn::Instance.

◆ storeWaterKin()

void moordyn::Line::storeWaterKin ( real  dt,
std::vector< std::vector< moordyn::real >>  zeta,
std::vector< std::vector< moordyn::real >>  f,
std::vector< std::vector< vec >>  u,
std::vector< std::vector< vec >>  ud 
)

store wave/current kinematics time series for this line

This is used when nodal approaches are selected, i.e. WaveKin = WAVES_FFT_NODE or WAVES_NODE, Currents = CURRENTS_STEADY_NODE or CURRENTS_DYNAMIC_NODE

Parameters
dtTime step
zetaThe wave elevations
fThe fluid fractions
uThe flow velocity
udThe flow acceleration
Exceptions
invalid_value_errorIf zeta, f, u and ud have not the same size, in both dimensions.
invalid_value_errorIf the length of zeta, f, u or ud is not equal to moordyn::Line::getN() + 1
Note
Working in progress

◆ updateUnstretchedLength()

void moordyn::Line::updateUnstretchedLength ( const moordyn::real  dt = 0.0)
inline

Update the unstretched length of the line, according to the velocity.

Parameters
dtThe time step. If zero, the initial unstretched length is set
Note
This function should be called after moordyn::Line::initialize()
Warning
The lines damping is not changed, which might affect the stability

Member Data Documentation

◆ lineId

size_t moordyn::Line::lineId

a line id which is guaranteed to be contiguous up from zero for the lines


The documentation for this class was generated from the following files: