Marine systems simulation
Ship::ShipSimObject Class Reference
+ Inheritance diagram for Ship::ShipSimObject:
+ Collaboration diagram for Ship::ShipSimObject:

Public Member Functions

 ShipSimObject (const string simObjectName, ISimObjectCreator *const creator)
 
void FinalSetup (const double T, const double *const X, ISimObjectCreator *const creator)
 
void GetCurrent (const double T, const double *const X, const double position[3], double *currentVelocity)
 
virtual void PreOdeFcn (const double T, const double *const X, IStateUpdater *updater)
 
virtual void OdeFcn (const double T, const double *const X, double *const XDot, const bool bIsMajorTimeStep)
 
const double * PositionNED (const double T, const double *const X)
 
const double * PositionZeroCrossNED (const double T, const double *const X)
 
const double * PositionAPNED (const double T, const double *const X)
 
const double * VelocityNED (const double T, const double *const X)
 
const double * QuaternionNED (const double T, const double *const X)
 
const double * EulerAnglesNED (const double T, const double *const X)
 
const double * RotationBody (const double T, const double *const X)
 
const double * OmegaNED (const double T, const double *const X)
 
void ComputeRollPitchYaw (const double T, const double *const X)
 
void ComputeRelativeVelocity (const double T, const double *const X)
 
const double * BodyForce (const double T, const double *const X)
 
const double * VelocityBody (const double T, const double *const X)
 
const double * FirstOrderWaveBody (const double T, const double *const X)
 The first order wave force, or froude-krylov force, as experienced in the body frame of the vessel The first order wave loads are defined as a sum of harmonic functions with the same frequency distribution as the wave elevation. More...
 
const double * WaveDriftBody (const double T, const double *const X)
 The mean value of the higher order wave forces as experienced in the body frame of the vessel. The mean value of the heave, roll and pitch direction equals zero. More...
 
const double *const FluidMemForceBody (const double T, const double *const X)
 
const double * RestoringForceBody (const double T, const double *const X)
 
const double * RelativeVelocityBody (const double T, const double *const X)
 
const double * ViscousForceBody (const double T, const double *const X)
 
const double * Encounter (const double T, const double *const X)
 
const double * CoriolisRigidBody (const double T, const double *const X)
 
const double * CoriolisFluidBody (const double T, const double *const X)
 
const double * SpeedInfo (const double T, const double *const X)
 
const double * RotationForceFrame (const double T, const double *const X)
 
const double * KineticEnergy (const double T, const double *const X)
 
const double * FluidEnergy (const double T, const double *const X)
 
const double * FluidEnergy2 (const double T, const double *const X)
 
const double * UFluidMemSpeed (const double T, const double *const X)
 
const double * FluidForceBySpeedBody (const double T, const double *const X, int index)
 
const double * FluidMemoryForceBody (const double T, const double *const X)
 
const double * CrossFlowBody (const double T, const double *const X)
 
const double * ResistanceBody (const double T, const double *const X)
 The calm water resistance of the hull in the surge direction The resistance force is the defined as the calm water resistance of the hull. The calm water resistance is divided into the Reynolds number ( \( R_n\)) dependent friction contribution, and the Froude number ( \( F_n\)) dependent wave generation from the flow pattern around the hull. The resistance is calculated by scaling the dimensional total resistance coefficient, \( C_{ts}\), by the wet surface, \( S\), on the current waterline. The total resistance coefficient is both \( R_n\) and \( F_n\) dependent and is defined as. More...
 
const double * AdditionalDampeningBody (const double T, const double *const X)
 
const double * ExternalForceNED (const double T, const double *const X)
 
const double * ExternalForceBody (const double T, const double *const X)
 
const double * ViscousRoll (const double T, const double *const X)
 
const double * SteadyStateControlBODY (const double T, const double *const X)
 
const double * SteadyStateControlNED (const double T, const double *const X)
 
const double * DriftSpeed (const double T, const double *const X)
 
- Public Member Functions inherited from Ship::ViscousShip
void SetupFromFile (std::string fname, EFileType type, HullData *H)
 
void SetHullData (HullData *HD)
 
void SetTheoryFormulation (int theory)
 
void SetPhiAmplitudeHorizon (double horizon)
 
void GetViscousAddedMass (double u, double v, double r, Eigen::Matrix3d &Mat)
 
double GetViscousSurge (double u, double v, double r)
 
double GetViscousSway (double u, double v, double r)
 
double GetViscousYaw (double u, double v, double r)
 
double GetViscousRoll (double phiDot, double u=0)
 
double EddieDamping (double phiDot, double u=0)
 
double BilgeKeelDamping (double phiDot, double u=0)
 
double FricitonDamping (double phiDot, double u=0)
 
double LiftDamping (double phiDot, double u=0)
 
double SkegDamping (double phiDot, double u=0)
 
double SkegLiftDamping (double phiDot, double u)
 
void AddToRollTimeHistory (const double T, double RollAngle, const bool bIsMajorTimeStep)
 
void MakeBildgeKeel (double m_BilgeKeelStartPosition, double m_BilgeKeelEndPosition, double m_BilgeKeelHeight)
 
void MakeSkeg ()
 
void SetHorizontalDampingBlend (double blendSpeed, double blendFraction=0.5)
 
void SetDpDampingX (double Xu, double Xv, double Xr)
 
void SetDpDampingY (double Yu, double Yv, double Yr)
 
void SetDpDampingN (double Nu, double Nv, double Nr)
 
double ManeuverReferenceSpeed () const
 
void ManeuverReferenceSpeed (double val)
 

Public Attributes

 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
 
double m_ViscousRoll [6]
 
double m_waveElevation
 
double m_waveSlope
 

Additional Inherited Members

- Public Types inherited from Ship::ViscousShip
enum  EFileType { FILE_TYPE_VERES , FILE_TYPE_WAMIT }
 
- Protected Member Functions inherited from Ship::ViscousShip
double deltaWingKP (double A_ratio)
 
double deltaWingKV (double A_ratio)
 
double deltaWingCL (double alfa, double Kp, double Kv)
 
void ReadMainDimensionsVeres (std::string fname)
 
void ReadMainDimensionsWamit (std::string fname)
 
- Protected Attributes inherited from Ship::ViscousShip
double m_maneuverReferenceSpeed
 
double m_viscousDampingBlender
 
double Yv_man
 
double Yr_man
 
double Yrd_man
 
double Yvd_man
 
double Nv_man
 
double Nr_man
 
double Nrd_man
 
double Nvd_man
 
double Xrr_man
 
double Xud_man
 
double Xvr_man
 
double Xvv_man
 
double Xvvvv_man
 
double Xu_dp
 
double Yu_dp
 
double Nu_dp
 
double Xv_dp
 
double Yv_dp
 
double Nv_dp
 
double Xr_dp
 
double Yr_dp
 
double Nr_dp
 
int m_bilgeKeelStartSection
 
int m_bilgeKeelEndSection
 
double m_startSectionFractionBK
 
double m_endSectionFractionBK
 
int rollDampingTheory
 
int m_withSkeg
 
int SkegStartSection
 
int SkegEndSection
 
int * SkegStartPoints
 
int * SkegEndPoints
 
HullDatamyHullData
 
double surfaceOfHull
 
double kinematicViscosity
 
double * SectionalNormalBilgeKeelForce
 
double * SectionalPressureOnHullFromBilgeKeels
 
double * SectionalEddieForce
 
double * bilgeKeelHeight
 
double * Ds
 
double * Bs
 
double * H0
 
double * SigmaS
 
double * A1
 
double * A3
 
double * Ms
 
double Lpp
 
double B
 
double T
 
double Cb
 
double VCG
 
double Lo
 
double Lr
 
double Kn
 
double * SectionLength
 
std::deque< double > phi_history
 
std::deque< double > time_history
 
double phiAmplitude
 
double phiAmplitudeTime
 
double m_rf
 
double OG
 
double phiAmplitudeHorizon
 
int numberOfStations
 
bool m_skegLift
 
bool m_hullLift
 

Member Function Documentation

◆ AdditionalDampeningBody()

const double * Ship::ShipSimObject::AdditionalDampeningBody ( const double  T,
const double *const  X 
)

The model supports additional user defined dampening forces. The forces are calculated from dimensionless dampening coefficients made dimensionless by \(m\) and by \(m\cdot L_{pp}\) for linear and rotational degrees of freedom respectively. The dampening force contains a linear and a quadratic modulus term which are controlled by separate dampening coefficients.

\[ F = C_{d_linear} \cdot m \cdot v + C_{d_quadratic} \cdot m \cdot v|v|\]

\[ M = C_{d_linear} \cdot m \cdot L_{pp} \cdot v + C_{d_quadratic} \cdot m \cdot L_{pp}\cdot v|v|\]

The dampening force decreases with increasing forward velocity by weighting the force with \(exp{-\frac{1}{2}|U|}\)

◆ BodyForce()

const double * Ship::ShipSimObject::BodyForce ( const double  T,
const double *const  X 
)

The total force acting on the vessel

◆ CoriolisFluidBody()

const double * Ship::ShipSimObject::CoriolisFluidBody ( const double  T,
const double *const  X 
)

Coriolis force of the added mass and inertia tensors as formulated in the center of gravity

◆ CoriolisRigidBody()

const double * Ship::ShipSimObject::CoriolisRigidBody ( const double  T,
const double *const  X 
)

Coriolis force of the rigid body mass and inertia tensors as formulated in the center of gravity

◆ CrossFlowBody()

const double * Ship::ShipSimObject::CrossFlowBody ( const double  T,
const double *const  X 
)

Cross flow are calculated by summation of section wise cross flow drag along the hull. Dampening in yaw is computed by calculating the moment contribution around the center of gravity for each sections cross flow drag force

\[ v_r = v-x\cdot r\]

\[ Y_{cross-flow} = -\frac{1}{2}\rho\int_{-Lpp/2}^{Lpp/2}T(x)C_d^2d(x)|v_r|v_r dx \]

\[ N_{cross-flow} = -\frac{1}{2}\rho\int_{-Lpp/2}^{Lpp/2}T(x)\cdot x\cdot C_d^2d(x)|v_r|v_r dx \]

◆ Encounter()

const double * Ship::ShipSimObject::Encounter ( const double  T,
const double *const  X 
)

Encounter frequency of the last wave in the sea environment list

◆ ExternalForceBody()

const double * Ship::ShipSimObject::ExternalForceBody ( const double  T,
const double *const  X 
)

External force in Body

◆ ExternalForceNED()

const double * Ship::ShipSimObject::ExternalForceNED ( const double  T,
const double *const  X 
)

External force in NED

◆ FirstOrderWaveBody()

const double * Ship::ShipSimObject::FirstOrderWaveBody ( const double  T,
const double *const  X 
)

\[ \tau_{1. order} = \sum_{n=0}^{n_{waves}}\zeta_n A(\omega_n,\psi_n,U) \cdot cos( \omega_{n} - k_n x + \epsilon_n) \]

◆ FluidEnergy()

const double * Ship::ShipSimObject::FluidEnergy ( const double  T,
const double *const  X 
)

The squared norm of the collection of fluid memory model states

◆ FluidMemForceBody()

const double *const Ship::ShipSimObject::FluidMemForceBody ( const double  T,
const double *const  X 
)

The dampening force from wave generation by the vessel motions. This is the output of the "fluid memory" model

◆ FluidMemoryForceBody()

const double * Ship::ShipSimObject::FluidMemoryForceBody ( const double  T,
const double *const  X 
)

Fluid memory forces on body

◆ KineticEnergy()

const double * Ship::ShipSimObject::KineticEnergy ( const double  T,
const double *const  X 
)

The kinetic energy of the rigid body of the vessel

◆ OmegaNED()

const double * Ship::ShipSimObject::OmegaNED ( const double  T,
const double *const  X 
)

The rotational speed of the ship in the global "north-east-down" coordinate system

◆ PositionAPNED()

const double * Ship::ShipSimObject::PositionAPNED ( const double  T,
const double *const  X 
)

The position of the intersection between the aft perpendicular, baseline in the global "north-east-down" coordinate system

◆ PositionNED()

const double * Ship::ShipSimObject::PositionNED ( const double  T,
const double *const  X 
)

The position of the center of gravity in the global "north-east-down" coordinate system

◆ PositionZeroCrossNED()

const double * Ship::ShipSimObject::PositionZeroCrossNED ( const double  T,
const double *const  X 
)

The position of the midpoint of the water plane at Lpp/2 in the global "north-east-down" coordinate system

◆ QuaternionNED()

const double * Ship::ShipSimObject::QuaternionNED ( const double  T,
const double *const  X 
)

The rotation quaternion of the body as w,x,y,z

◆ RelativeVelocityBody()

const double * Ship::ShipSimObject::RelativeVelocityBody ( const double  T,
const double *const  X 
)

The relative velocity of the ship relative to the ambient (current) velocity

◆ ResistanceBody()

const double * Ship::ShipSimObject::ResistanceBody ( const double  T,
const double *const  X 
)

\[ C_{ts} = (1+k) (C_f(R_n) + \Delta C_f(R_n)) + C_r(F_n) \]

Where

  • \(C_{ts}\) is the total resistance coefficient
  • \(k\) is a form factor modification of the friction resistance due to the difference between the a flat plate and a hull form defined as k= \((0.6+145(\frac{C_b}{L_{wl}}\sqrt{B\cdot(T_a+T_f)})^(2.5))\frac{C_b}{L_{wl}}\sqrt{B\cdot(T_a+T_f)}\)
    • \({C_b}\) is the block coefficient
    • \(B\) is the breadth
    • \(T_a\) and \(T_f\) is the aft and fore draught
  • \(C_f\) is the \(R_n\) dependent friction coefficient of a flat plate defined as \(C_f = \frac{0.0075}{(log_{10}R_n-2)^2}\)
  • \(\Delta C_f\) is the additional friction resistance due to roughness of the hull
  • \(C_r\) is the \(F_n\) dependent resistance from generation of waves The resistance due to wind is neglected in this resistance calculation to separate the wind force from the resistance. The total resistance is calculated from

    \[ R_{ts} = \frac{1}{2}\rho U^2 \cdot S \cdot C_{ts}\]

  • \(\rho\) is the density of sea water
  • \(U\) is the surge speed of the vessel

◆ RestoringForceBody()

const double * Ship::ShipSimObject::RestoringForceBody ( const double  T,
const double *const  X 
)

The hydrostatic restoring force on the body

◆ RotationBody()

const double * Ship::ShipSimObject::RotationBody ( const double  T,
const double *const  X 
)

The vessels roll, pitch and yaw angles in the local coordinate system. The yaw angle is relative to the north direction.

◆ RotationForceFrame()

const double * Ship::ShipSimObject::RotationForceFrame ( const double  T,
const double *const  X 
)

Rotation body if constrained to move parallel to the still ocean surface

◆ SpeedInfo()

const double * Ship::ShipSimObject::SpeedInfo ( const double  T,
const double *const  X 
)

Information about the reference speeds in the potential data and the current speed of the vessel

◆ SteadyStateControlBODY()

const double * Ship::ShipSimObject::SteadyStateControlBODY ( const double  T,
const double *const  X 
)

Steady state control in BODY coordinates

◆ SteadyStateControlNED()

const double * Ship::ShipSimObject::SteadyStateControlNED ( const double  T,
const double *const  X 
)

Steady state control in NED coordinates

◆ VelocityBody()

const double * Ship::ShipSimObject::VelocityBody ( const double  T,
const double *const  X 
)

The full velocity vector of the vessel in body coordinates

◆ VelocityNED()

const double * Ship::ShipSimObject::VelocityNED ( const double  T,
const double *const  X 
)

The six-dof velocity vector for the center of gravity in the global "north-east-down" coordinate system

◆ ViscousForceBody()

const double * Ship::ShipSimObject::ViscousForceBody ( const double  T,
const double *const  X 
)

Viscous forces on body

◆ WaveDriftBody()

const double * Ship::ShipSimObject::WaveDriftBody ( const double  T,
const double *const  X 
)

\[ \tau_{wave-drift} = \sum_{n=0}^{n_{waves}} {\zeta_{n}}^2 A(\omega_n,\psi_n,U) \]


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