beams Namespace Reference

Kynema API: kynema::beams Namespace Reference
Kynema API
A flexible multibody structural dynamics code for wind turbines
Loading...
Searching...
No Matches
kynema::beams Namespace Reference

Classes

struct  CalculateForceFC
 
struct  CalculateForceFD
 
struct  CalculateInertialQuadraturePointValues
 
struct  CalculateJacobian
 Functor to calculate Jacobians and unit tangent vectors at quadrature points for beam elements. More...
 
struct  CalculateOuu
 
struct  CalculatePuu
 
struct  CalculateQPPosition
 Functor to calculate current position and orientation at quadrature points. More...
 
struct  CalculateQuadraturePointValues
 
struct  CalculateQuu
 
struct  CalculateStiffnessQuadraturePointValues
 
struct  CalculateStrain
 
struct  CalculateSystemMatrix
 
struct  CalculateTemporaryVariables
 
struct  HollowCircleProperties
 Struct containing geometric properties for a hollow circular cross-section. More...
 
struct  IntegrateInertiaMatrixElement
 
struct  IntegrateResidualVectorElement
 
struct  IntegrateStiffnessMatrixElement
 
struct  InterpolateQPPosition
 Interpolates quadrature point positions from nodal positions using shape functions. More...
 
struct  InterpolateQPRotation
 A Kernel which interpolates a rotation quaternion on a given element from its nodes to all of it quadrature points. More...
 
struct  InterpolateQPState_r
 Interpolates the rotation (r) part of the state at a quadrature point. More...
 
struct  InterpolateQPState_rprime
 Interpolates the rotation derivative (r') part of the state at a quadrature point. More...
 
struct  InterpolateQPState_u
 Interpolates the displacement (u) part of the state at a quadrature point. More...
 
struct  InterpolateQPState_uprime
 Interpolates the displacement derivative (u') part of the state at a quadrature point. More...
 
struct  InterpolateQPVector
 A Kernel which interpolates a vector quantity from nodes on a given element to a quadrature point given a basis function. More...
 
struct  InterpolateToQuadraturePointForInertia
 
struct  InterpolateToQuadraturePointForStiffness
 
struct  InterpolateToQuadraturePoints
 Interpolates various quantities from nodes to quadrature points for beam elements. More...
 
struct  UpdateNodeState
 
struct  UpdateNodeStateElement
 

Functions

std::vector< std::array< double, 2 > > CreateTrapezoidalQuadrature (std::span< const double > grid)
 Creates a trapezoidal quadrature rule based on a given grid.
 
std::vector< std::array< double, 2 > > CreateGaussLegendreLobattoQuadrature (std::span< const double > grid, size_t order)
 
std::vector< std::array< double, 2 > > CreateGaussLegendreQuadrature (const size_t order)
 Creates Gauss-Legendre (GL) quadrature points and weights on [-1, 1].
 
static std::array< std::array< double, 6 >, 6 > GenerateStiffnessMatrix (double EA, double EI_x, double EI_y, double GKt, double GA, double kxs, double kys, double x_C, double y_C, double theta_p, double x_S, double y_S, double theta_s)
 Generates a 6x6 cross-sectional stiffness matrix for use in beam elements.
 
static std::array< std::array< double, 6 >, 6 > GenerateMassMatrix (double m, double I_x, double I_y, double I_p, double x_G=0., double y_G=0., double theta_i=0.)
 Generates a 6x6 cross-sectional mass matrix for use in beam elements.
 
static HollowCircleProperties CalculateHollowCircleProperties (double outer_diameter, double wall_thickness, double nu=0.33)
 Calculates geometric properties for a hollow circular cross-section.
 
static BeamSection GenerateHollowCircleSection (double s, double E, double G, double rho, double outer_diameter, double wall_thickness, double nu, double x_C=0., double y_C=0., double theta_p=0., double x_S=0., double y_S=0., double theta_s=0., double x_G=0., double y_G=0., double theta_i=0)
 Generates a BeamSection with 6x6 mass and stiffness matrices for a hollow circular cross-section.
 
template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void InterpVector3 (const shape_matrix_type &shape_matrix, const node_type &node_v, const qp_type &qp_v)
 
template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void InterpVector4 (const shape_matrix_type &shape_matrix, const node_type &node_v, const qp_type &qp_v)
 
template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void InterpQuaternion (const shape_matrix_type &shape_matrix, const node_type &node_v, const qp_type &qp_v)
 
template<typename shape_matrix_type , typename jacobian_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void InterpVector3Deriv (const shape_matrix_type &shape_matrix_deriv, const jacobian_type &jacobian, const node_type &node_v, const qp_type &qp_v)
 
template<typename shape_matrix_type , typename jacobian_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void InterpVector4Deriv (const shape_matrix_type &shape_matrix_deriv, const jacobian_type &jacobian, const node_type &node_v, const qp_type &qp_v)
 
void PopulateNodeX0 (const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double *[7], Kokkos::LayoutStride, Kokkos::HostSpace > &node_x0)
 Populate the node initial position and orientation.
 
void PopulateQPWeight (const BeamElement &elem, const Kokkos::View< double *, Kokkos::LayoutStride, Kokkos::HostSpace > &qp_weight)
 Populate the integration weights at each quadrature point.
 
std::vector< double > MapNodePositions (const BeamElement &elem, std::span< const Node > nodes)
 Map node positions from [0,1] to [-1,1].
 
void PopulateShapeFunctionValues (const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &shape_interp)
 Populate shape function values at each quadrature point.
 
void PopulateShapeFunctionDerivatives (const BeamElement &elem, std::span< const Node > nodes, const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &shape_deriv)
 Populate shape function derivatives at each quadrature point.
 
std::vector< double > MapSectionPositions (const BeamElement &elem)
 Map section positions from [0,1] to [-1,1].
 
void PopulateQPMstar (const BeamElement &elem, const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &qp_Mstar)
 Populate mass matrix values at each quadrature point.
 
void PopulateQPCstar (const BeamElement &elem, const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &qp_Cstar)
 Populate stiffness matrix values at each quadrature point.
 

Function Documentation

◆ CalculateHollowCircleProperties()

static HollowCircleProperties kynema::beams::CalculateHollowCircleProperties ( double  outer_diameter,
double  wall_thickness,
double  nu = 0.33 
)
static

Calculates geometric properties for a hollow circular cross-section.

Parameters
outer_diameterOuter diameter of the hollow circle [Length]
wall_thicknessWall thickness [Length]
nuPoisson's ratio for shear correction factor calculation [dimensionless]
Returns
HollowCircleProperties struct containing geometric properties

◆ CreateGaussLegendreLobattoQuadrature()

std::vector< std::array< double, 2 > > kynema::beams::CreateGaussLegendreLobattoQuadrature ( std::span< const double >  grid,
size_t  order 
)
inline

◆ CreateGaussLegendreQuadrature()

std::vector< std::array< double, 2 > > kynema::beams::CreateGaussLegendreQuadrature ( const size_t  order)
inline

Creates Gauss-Legendre (GL) quadrature points and weights on [-1, 1].

GL quadrature provides optimal accuracy for polynomial integration. The points are roots of P_n(x) where P_n is the nth Legendre polynomial. GL quadrature does NOT include the endpoints (-1, 1).

Parameters
orderNumber of quadrature points (n >= 1). Returns n quadrature points.
Returns
Vector of {point, weight} pairs for GL quadrature, sorted by point value

◆ CreateTrapezoidalQuadrature()

std::vector< std::array< double, 2 > > kynema::beams::CreateTrapezoidalQuadrature ( std::span< const double >  grid)
inline

Creates a trapezoidal quadrature rule based on a given grid.

This function generates a set of quadrature points and weights for trapezoidal integration over a specified grid. The quadrature points are mapped to the interval [-1, 1].

Parameters
gridA span of grid points defining the integration domain.
Returns
A vector of arrays, where each array contains a quadrature point and its corresponding weight.

◆ GenerateHollowCircleSection()

static BeamSection kynema::beams::GenerateHollowCircleSection ( double  s,
double  E,
double  G,
double  rho,
double  outer_diameter,
double  wall_thickness,
double  nu,
double  x_C = 0.,
double  y_C = 0.,
double  theta_p = 0.,
double  x_S = 0.,
double  y_S = 0.,
double  theta_s = 0.,
double  x_G = 0.,
double  y_G = 0.,
double  theta_i = 0 
)
static

Generates a BeamSection with 6x6 mass and stiffness matrices for a hollow circular cross-section.

This is a convenience function specifically for hollow circular sections commonly used in wind turbine towers. It calculates the geometric properties and then calls the general GenerateMassMatrix and GenerateStiffnessMatrix functions.

Parameters
sNormalized position along beam element [dimensionless]
EYoung's modulus [Force/Length²]
GShear modulus [Force/Length²]
rhoMaterial density [Mass/Length³]
outer_diameterOuter diameter of the hollow circle [Length]
wall_thicknessWall thickness [Length]
nuPoisson's ratio [dimensionless]
x_Cx-coordinate of elastic centroid relative to reference point [Length]
y_Cy-coordinate of elastic centroid relative to reference point [Length]
theta_pRotation angle from reference axes to principal bending axes [radians]
x_Sx-coordinate of shear center relative to reference point [Length]
y_Sy-coordinate of shear center relative to reference point [Length]
theta_sRotation angle from reference axes to principal shear axes [radians]
x_Gx-coordinate of center of gravity relative to reference point [Length]
y_Gy-coordinate of center of gravity relative to reference point [Length]
theta_iRotation angle from reference axes to principal inertia axes [radians]
Returns
BeamSection object containing stiffness matrix, mass matrix, and position
Note
For hollow circular sections, the elastic centroid, shear center, and center of gravity all coincide at the geometric center, so offset parameters are typically zero.
Principal axes align with the reference axes due to circular symmetry.

◆ GenerateMassMatrix()

static std::array< std::array< double, 6 >, 6 > kynema::beams::GenerateMassMatrix ( double  m,
double  I_x,
double  I_y,
double  I_p,
double  x_G = 0.,
double  y_G = 0.,
double  theta_i = 0. 
)
static

Generates a 6x6 cross-sectional mass matrix for use in beam elements.

This function constructs the cross-sectional mass matrix that relates the generalized accelerations to the generalized inertial forces in a beam cross-section. Returns the mass matrix at a given point 'O' and w.r.t. given orientation axes based on the values at the center of gravity 'G' and in the inertial axis frame.

Parameters
mMass per unit length [Mass/Length]
I_xMass moment of inertia about local x-axis in principal/inertial frame and at center of gravity G [Mass×Length]
I_yMass moment of inertia about local y-axis in principal/inertial frame and at center of gravity G [Mass×Length]
I_pPolar mass moment of inertia in principal/inertial frame and at center of gravity G (I_x + I_y) [Mass×Length]
x_Gx-coordinate of center of gravity relative to reference point i.e. distance FROM O TO G [Length]
y_Gy-coordinate of center of gravity relative to reference point i.e. distance FROM O TO G [Length]
theta_iRotation angle (around z) FROM reference axes TO principal inertial axes [radians]
Returns
6x6 cross-sectional mass matrix

◆ GenerateStiffnessMatrix()

static std::array< std::array< double, 6 >, 6 > kynema::beams::GenerateStiffnessMatrix ( double  EA,
double  EI_x,
double  EI_y,
double  GKt,
double  GA,
double  kxs,
double  kys,
double  x_C,
double  y_C,
double  theta_p,
double  x_S,
double  y_S,
double  theta_s 
)
static

Generates a 6x6 cross-sectional stiffness matrix for use in beam elements.

This function constructs the cross-sectional stiffness matrix that relates the generalized strains to the generalized forces/moments in a beam cross-section. The matrix accounts for coupling between axial, bending, shear, and torsional deformations due to offset between the elastic centroid and shear center locations. Returns a 6x6 stiffness matrix at the cross section origin 'O', based on inputs at the centroid 'C' and shear center 'S' locations.

Parameters
EAAxial stiffness (E×A) [Force]
EI_xBending stiffness around local x-axis in principal frame and at the centroid 'C' [Force × Length²]
EI_yBending stiffness around local y-axis in principal frame and at the centroid 'C' [Force × Length²]
GKtTorsional stiffness around principal frame and at shear center 'S' (G×J) [Force × Length²]
GAShear stiffness around principal frame and at shear center 'S' (G×A) [Force]
kxsShear correction factor in x-direction [dimensionless]
kysShear correction factor in y-direction [dimensionless]
x_Cx-coordinate of elastic centroid relative to reference point i.e. distance FROM 'O' TO 'C' [Length]
y_Cy-coordinate of elastic centroid relative to reference point i.e. distance FROM 'O' TO 'C' [Length]
theta_pRotation angle (around z) FROM reference axes TO principal bending axes [radians]
x_Sx-coordinate of shear center relative to reference point i.e. distance FROM 'O' TO 'S' [Length]
y_Sy-coordinate of shear center relative to reference point i.e. distance FROM 'O' TO 'S' [Length]
theta_sRotation angle (around z) FROM reference axes TO principal shear axes [radians]
Returns
6x6 cross-sectional stiffness matrix

◆ InterpQuaternion()

template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void kynema::beams::InterpQuaternion ( const shape_matrix_type &  shape_matrix,
const node_type &  node_v,
const qp_type &  qp_v 
)

◆ InterpVector3()

template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void kynema::beams::InterpVector3 ( const shape_matrix_type &  shape_matrix,
const node_type &  node_v,
const qp_type &  qp_v 
)

◆ InterpVector3Deriv()

template<typename shape_matrix_type , typename jacobian_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void kynema::beams::InterpVector3Deriv ( const shape_matrix_type &  shape_matrix_deriv,
const jacobian_type &  jacobian,
const node_type &  node_v,
const qp_type &  qp_v 
)

◆ InterpVector4()

template<typename shape_matrix_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void kynema::beams::InterpVector4 ( const shape_matrix_type &  shape_matrix,
const node_type &  node_v,
const qp_type &  qp_v 
)

◆ InterpVector4Deriv()

template<typename shape_matrix_type , typename jacobian_type , typename node_type , typename qp_type >
KOKKOS_INLINE_FUNCTION void kynema::beams::InterpVector4Deriv ( const shape_matrix_type &  shape_matrix_deriv,
const jacobian_type &  jacobian,
const node_type &  node_v,
const qp_type &  qp_v 
)

◆ MapNodePositions()

std::vector< double > kynema::beams::MapNodePositions ( const BeamElement elem,
std::span< const Node nodes 
)
inline

Map node positions from [0,1] to [-1,1].

◆ MapSectionPositions()

std::vector< double > kynema::beams::MapSectionPositions ( const BeamElement elem)
inline

Map section positions from [0,1] to [-1,1].

◆ PopulateNodeX0()

void kynema::beams::PopulateNodeX0 ( const BeamElement elem,
std::span< const Node nodes,
const Kokkos::View< double *[7], Kokkos::LayoutStride, Kokkos::HostSpace > &  node_x0 
)
inline

Populate the node initial position and orientation.

◆ PopulateQPCstar()

void kynema::beams::PopulateQPCstar ( const BeamElement elem,
const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &  qp_Cstar 
)
inline

Populate stiffness matrix values at each quadrature point.

◆ PopulateQPMstar()

void kynema::beams::PopulateQPMstar ( const BeamElement elem,
const Kokkos::View< double *[6][6], Kokkos::LayoutStride, Kokkos::HostSpace > &  qp_Mstar 
)
inline

Populate mass matrix values at each quadrature point.

◆ PopulateQPWeight()

void kynema::beams::PopulateQPWeight ( const BeamElement elem,
const Kokkos::View< double *, Kokkos::LayoutStride, Kokkos::HostSpace > &  qp_weight 
)
inline

Populate the integration weights at each quadrature point.

◆ PopulateShapeFunctionDerivatives()

void kynema::beams::PopulateShapeFunctionDerivatives ( const BeamElement elem,
std::span< const Node nodes,
const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &  shape_deriv 
)
inline

Populate shape function derivatives at each quadrature point.

◆ PopulateShapeFunctionValues()

void kynema::beams::PopulateShapeFunctionValues ( const BeamElement elem,
std::span< const Node nodes,
const Kokkos::View< double **, Kokkos::LayoutStride, Kokkos::HostSpace > &  shape_interp 
)
inline

Populate shape function values at each quadrature point.