|
MoorDyn
|
A mooring line. More...
#include <Line.hpp>


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 vec & | getNodePos (unsigned int i) const |
| Get the position of a node. More... | |
| const vec & | getNodeVel (unsigned int i) const |
| Get the velocity of a node. More... | |
| const vec & | getNodeForce (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 vec & | getNodeBendStiff (unsigned int i) const |
| Get the tension on a node, including the internal line damping. More... | |
| const vec & | getNodeWeight (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 real & | getNodeCurv (unsigned int i) const |
| Get the line curvature at a node position. More... | |
| const mat & | getNodeM (unsigned int i) const |
| Get the mass and added mass matrix. More... | |
| std::vector< vec > | getNodeCoordinates () 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, real > | getDrag () 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... | |
| Log * | GetLogger () 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. | |
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,
, the number of segments,
, and the material density and Young's modulus,
&
, a natural oscillation period can be defined:
![]()
The integration time step (moordyn::MoorDyn.dtM0) should be smaller than this natural period to avoid numerical instabilities
| moordyn::Line::Line | ( | moordyn::Log * | log, |
| size_t | lineId | ||
| ) |
Constructor.
| log | Logging handler |
| lineId | the incremental Id of the line |
| 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.
| firstNodeIdx | Index of the first node of the segment |
| secondNodeIdx | Index of the second node of the segment |
| surfaceHeight | Height of the water surface (assumed locally horizontal) |


|
virtual |
Unpack the data to restore the Serialized information.
This is the inverse of Serialize(void)
| data | The packed data |
Implements moordyn::io::IO.

|
inline |
Check if pressure bending forces are considered.
|
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.
|
inline |
Get the constant stiffness of the line.
This value is useless if non-linear stiffness is considered
Get the drag coefficients.
|
inline |
Get the elastic model used by the line.
Get the bending moment at the end point.
The bending moment is defined as
, with
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
| end_point | Either ENDPOINT_B or ENDPOINT_A |
| rod_end_point | Either ENDPOINT_B or ENDPOINT_A |
| invalid_value_error | If either end_point or rod_end_point are not valid end point qualifiers |

|
inline |
Get the force, moment and mass at the line endpoint.
| Fnet_out | The output force vector |
| Moment_out | The output moment vector |
| M_out | The output mass matrix |
| end_point | Either ENDPOINT_TOP or ENDPOINT_BOTTOM |
| invalid_value_error | If end_point is not a valid end point qualifier |
|
inline |
Get the tensions at the fairlead and anchor in a FASTv7 friendly way.
| FairHTen | Horizontal tension on the fairlead |
| FairVTen | Vertical tension on the fairlead |
| AnchHTen | Horizontal tension on the anchor |
| AnchVTen | Vertical tension on the anchor |
| real moordyn::Line::GetLineOutput | ( | OutChanProps | outChan | ) |
Get line output.
This function is useful when outputs are set in the line properties
| outChan | The output channel/field |

|
inline |
Number of segments.
The number of nodes can be computed as moordyn::Line::getN() + 1

|
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
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the array of coordinates of all nodes along the line.
|
inline |
Get the line curvature at a node position.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the drag force acting on the node.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
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()
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |

|
inline |
Get the Froude-Krilov force acting on the node.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the mass and added mass matrix.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the position of a node.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the sea bed reaction force acting on the node.
This is eventually including the friction force
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
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.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |

|
inline |
Get the velocity of a node.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
inline |
Get the weight and bouyancy force acting on the node.
| i | The line node index |
| invalid_value_error | If the node index i is bigger than the number of nodes, moordyn::Line::N + 1 |
|
virtual |
Calculate forces and get the derivative of the line's states.
| drdt | Output state derivatives |
| nan_error | If nan values are detected in any node position |
Implements moordyn::Instance.

|
inline |
Get the unstretched length of the line.
|
inline |
Get the VTK writer.
This function is useful for writing multiblock .vtm files
| void moordyn::Line::initialize | ( | ) |
Initialize the line object.
| moordyn::output_file_error | If an outfile has been provided, but it cannot be written |
| invalid_value_error | If there is no enough water depth |

|
inlinevirtual |
Sets the initial line state.
| state | The line state for initializing (see Line::setState for the state structure0) |
Implements moordyn::Instance.
|
inline |
Get whether the line is governed by a non-linear stiffness or a constant one.
| void moordyn::Line::saveVTK | ( | const char * | filename | ) |
Save the line on a VTK (.vtp) file.
| filename | The output file name |
| output_file_error | If If the file cannot be saved |

|
inline |
Scale the drag coefficients.
| scaler | The drag coefficients scale factor |
|
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
Implements moordyn::io::IO.

|
inline |
Set the constant stiffness of the line.
This value is useless if non-linear stiffness is considered
| EA | The constant axial stiffness EA value (N) |
Set the drag coefficients.
| cdn | Normal (transversal) coefficient |
| cdt | tangential (axial) coefficient |
Set the position and velocity of an end point.
| r | Position |
| rd | Velocity |
| end_point | Either ENDPOINT_TOP or ENDPOINT_BOTTOM |
| invalid_value_error | If end_point is not a valid end point qualifier |
set end node unit vector
This method is called by an eventually attached Rod, only applicable for bending stiffness
| q | The direction unit vector |
| end_point | Either ENDPOINT_B or ENDPOINT_A |
| rod_end_point | Either ENDPOINT_B or ENDPOINT_A |
| invalid_value_error | If either end_point or rod_end_point are not valid end point qualifiers |
|
inline |
Set the initial Tmax values.
| Tmax0 | Initial Tmax value |
|
inline |
Set the initial mean tension values.
| Tmean0 | Initial Tmean value |
| 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.
| p | The pressure values (Pa) at each node |
| invalid_value_error | If p has wrong size |
|
virtual |
Set the line state.
| r | The line state matrix. See Line::setState in Line.cpp for structure |
Implements moordyn::Instance.
|
inline |
Set the line simulation time.
| time | Simulation time |
|
inline |
Set the unstretched length of the line.
This function is consistently changing the unstretched length of each segment, moordyn::Line::l
| len | The unstretched length, moordyn::Line::UnstrLen |
|
inline |
Set the unstretched length rate of change of the line.
| v | The unstretched length rate of change, moordyn::Line::UnstrLend |
| 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.
| number | Line number |
| props | Line properties |
| l | Unstretched line length |
| n | Number of segments |
| env_in | Global struct that holds environmental settings |
| outfile | The outfile where information shall be written |
| channels | The channels/fields that shall be printed in the file |
| dtM0 | moordyn internal timestep. Used in VIV |

|
inline |
|
inlinevirtual |
Get the dimension of the state variable.
Reimplemented from moordyn::Instance.
|
inlinevirtual |
Get the number of state variables required by this instance.
Reimplemented from moordyn::Instance.
| 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
| dt | Time step |
| zeta | The wave elevations |
| f | The fluid fractions |
| u | The flow velocity |
| ud | The flow acceleration |
| invalid_value_error | If zeta, f, u and ud have not the same size, in both dimensions. |
| invalid_value_error | If the length of zeta, f, u or ud is not equal to moordyn::Line::getN() + 1 |
|
inline |
Update the unstretched length of the line, according to the velocity.
| dt | The time step. If zero, the initial unstretched length is set |
| size_t moordyn::Line::lineId |
a line id which is guaranteed to be contiguous up from zero for the lines