Private API

Math

GXBeam.rotation_parameter_scalingFunction
rotation_parameter_scaling(θ)

Extract a scaling parameter which may be multiplied by the angular parameters to yield the Wiener-Milenkovic rotation parameters. Use of this scaling parameter allows deflections greater than 360 degrees.

source
GXBeam.get_CFunction
get_C(θ)

Returns the transformation matrix C given the three angular parameters in θ.

source
GXBeam.get_C_θFunction
get_C_θ([C, ] θ)

Calculate the derivative of the Wiener-Milenkovic transformation matrix C with respect to each of the rotation parameters in θ.

source
GXBeam.get_QFunction
get_Q(θ)

Calculate the matrix Q as defined in the paper "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu given the rotational parameters in θ.

source
GXBeam.get_Q_θFunction
get_Q_θ(θ)
get_Q_θ(Q, θ)

Calculate the derivative of the matrix Q with respect to each of the rotation parameters in θ.

source
GXBeam.get_QinvFunction
get_Qinv(θ)

Calculate the matrix inverse Qinv as defined in the paper "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu given the rotational parameters in θ.

source
GXBeam.get_Qinv_θFunction
get_Qinv_θ(θ)

Calculate the derivative of the matrix inverse Qinv with respect to each of the rotation parameters in θ.

source
GXBeam.get_ΔQFunction
get_ΔQ(θ, Δθ [, Q])

Calculate the matrix ΔQ for structural damping calculations

source
GXBeam.get_ΔQ_θFunction
get_ΔQ_θ(θ, Δθ, [Q, Q_θ1, Q_θ2, Q_θ3])

Calculate the derivative of the matrix ΔQ with respect to each of the rotation parameters in θ.

source
GXBeam.mul3Function
mul3(A_1, A_2, A_3, b)

Return the product of a 3x3x3 tensor represented by A_1, A_2, and A_3 with the vector b.

source

Body Frame

GXBeam.update_body_acceleration_indices!Function
update_body_acceleration_indices!(system, prescribed_conditions)
update_body_acceleration_indices!(indices, prescribed_conditions)

Updates the state variable indices corresponding to the body frame accelerations to correspond to the provided prescribed conditions.

source
GXBeam.body_accelerationsFunction
body_accelerations(system, x=system.x; linear_acceleration=zeros(3), angular_acceleration=zeros(3))

Extract the linear and angular acceleration of the body frame from the state vector, if applicable. Otherwise return the provided linear and angular acceleration. This function is applicable only for steady state and initial condition analyses.

source

Points

GXBeam.point_loadsFunction
point_loads(x, ipoint, icol, force_scaling, prescribed_conditions)

Extract the loads F and M of point ipoint from the state variable vector or prescribed conditions.

source
GXBeam.point_load_jacobiansFunction
point_load_jacobians(x, ipoint, icol, force_scaling, prescribed_conditions)

Calculate the load jacobians F_θ, F_F, M_θ, and M_M of point ipoint.

source
GXBeam.point_displacementFunction
point_displacement(x, ipoint, icol_point, prescribed_conditions)

Extract the displacements u and θ of point ipoint from the state variable vector or prescribed conditions.

source
GXBeam.point_displacement_ratesFunction
point_displacement_rates(dx, ipoint, icol, prescribed_conditions)

Extract the displacement rates udot and θdot of point ipoint from the rate variable vector.

source
GXBeam.point_velocity_ratesFunction
point_velocity_rates(dx, ipoint, icol_point)

Extract the velocity rates Vdot and Ωdot of point ipoint from the rate variable vector.

source
GXBeam.point_velocitiesFunction
point_velocities(x, ipoint, icol_point)

Extract the velocities V and Ω of point ipoint from the state variable vector

source
GXBeam.initial_point_displacementFunction
initial_point_displacement(x, ipoint, icol_point, prescribed_conditions,
    rate_vars)

Extract the displacements u and θ of point ipoint from the state variable vector or prescribed conditions for an initial condition analysis.

source
GXBeam.initial_point_velocity_ratesFunction
initial_point_velocity_rates(x, ipoint, icol_point, prescribed_conditions,
    Vdot0, Ωdot0, rate_vars)

Extract the velocity rates Vdot and Ωdot of point ipoint from the state variable vector or provided initial conditions. Note that Vdot and Ωdot in this case do not include any contributions resulting from body frame motion.

source
GXBeam.initial_point_displacement_jacobianFunction
initial_point_displacement_jacobian(ipoint, icol_point, prescribed_conditions,
    rate_vars)

Extract the displacement jacobians u_u and θ_θ of point ipoint from the state variable vector or prescribed conditions for an initial condition analysis.

source
GXBeam.initial_point_velocity_rate_jacobianFunction
initial_point_velocity_rate_jacobian(ipoint, icol_point, prescribed_conditions,
    rate_vars)

Return the velocity rate jacobians Vdot_Vdot and Ωdot_Ωdot of point ipoint. Note that Vdot and Ωdot in this case do not include any contributions resulting from body frame motion.

source
GXBeam.static_point_propertiesFunction
static_point_properties(x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity)

Calculate/extract the point properties needed to construct the residual for a static analysis

source
GXBeam.steady_point_propertiesFunction
steady_point_properties(x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))

Calculate/extract the point properties needed to construct the residual for a steady state analysis

source
GXBeam.initial_point_propertiesFunction
initial_point_properties(x, indices, rate_vars, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
    u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate/extract the point properties needed to construct the residual for a time domain analysis initialization.

source
GXBeam.newmark_point_propertiesFunction
newmark_point_properties(x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    udot_init, θdot_init, Vdot_init, Ωdot_init, dt)

Calculate/extract the point properties needed to construct the residual for a newmark-scheme time stepping analysis

source
GXBeam.dynamic_point_propertiesFunction
dynamic_point_properties(dx, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)

Calculate/extract the point properties needed to construct the residual for a dynamic analysis

source
GXBeam.expanded_steady_point_propertiesFunction
expanded_steady_point_properties(x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))

Calculate/extract the point properties needed to construct the residual for a constant mass matrix system.

source
GXBeam.expanded_dynamic_point_propertiesFunction
expanded_dynamic_point_properties(dx, x, indices, force_scaling, assembly,
    ipoint, prescribed_conditions, point_masses, gravity, linear_velocity,
    angular_velocity)

Calculate/extract the point properties needed to construct the residual for a constant mass matrix system.

source
GXBeam.static_point_residual!Function
static_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity)

Calculate and insert the residual entries corresponding to a point for a static analysis into the system residual vector.

source
GXBeam.steady_point_residual!Function
steady_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    linear_acceleration, angular_acceleration)

Calculate and insert the residual entries corresponding to a point for a steady state analysis into the system residual vector.

source
GXBeam.initial_point_residual!Function
initial_point_residual!(resid, x, indices, rate_vars,
    force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
    u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate and insert the residual entries corresponding to a point for the initialization of a time domain analysis into the system residual vector.

source
GXBeam.newmark_point_residual!Function
newmark_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    udot_init, θdot_init, Vdot_init, Ωdot_init, dt)

Calculate and insert the residual entries corresponding to a point for a newmark-scheme time marching analysis into the system residual vector.

source
GXBeam.dynamic_point_residual!Function
dynamic_point_residual!(resid, dx, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)

Calculate and insert the residual entries corresponding to a point for a dynamic analysis into the system residual vector.

source
GXBeam.expanded_steady_point_residual!Function
expanded_steady_point_residual!(resid, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    linear_acceleration, angular_acceleration)

Calculate and insert the residual entries corresponding to a point into the system residual vector for a constant mass matrix system.

source
GXBeam.expanded_dynamic_point_residual!Function
expanded_dynamic_point_residual!(resid, dx, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)

Calculate and insert the residual entries corresponding to a point into the system residual vector for a constant mass matrix system.

source
GXBeam.static_point_jacobian_propertiesFunction
static_point_jacobian_properties(properties, x, indices, force_scaling, assembly,
    ipoint, prescribed_conditions, point_masses, gravity)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a static analysis

source
GXBeam.steady_point_jacobian_propertiesFunction
steady_point_jacobian_properties(properties, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity,
    ub_ub, θb_θb, vb_vb, ωb_ωb, ab_ab, αb_αb)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a steady state analysis

source
GXBeam.initial_point_jacobian_propertiesFunction
initial_point_jacobian_properties(properties, x, indices, rate_vars,
    force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
    u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a Newmark scheme time marching analysis

source
GXBeam.newmark_point_jacobian_propertiesFunction
newmark_point_jacobian_properties(properties, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity,
    udot_init, θdot_init, Vdot_init, Ωdot_init, dt)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a Newmark scheme time marching analysis

source
GXBeam.dynamic_point_jacobian_propertiesFunction
dynamic_point_jacobian_properties(properties, dx, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a dynamic analysis

source
GXBeam.expanded_steady_point_jacobian_propertiesFunction
expanded_steady_point_jacobian_properties(properties, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a constant mass matrix system

source
GXBeam.expanded_dynamic_point_jacobian_propertiesFunction
expanded_dynamic_point_jacobian_properties(properties, dx, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity)

Calculate/extract the point properties needed to calculate the jacobian entries corresponding to a point for a constant mass matrix system

source
GXBeam.mass_matrix_point_jacobian_propertiesFunction
mass_matrix_point_jacobian_properties(x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses)

Calculate/extract the point properties needed to calculate the mass matrix jacobian entries corresponding to a point

source
GXBeam.expanded_mass_matrix_point_jacobian_propertiesFunction
expanded_mass_matrix_point_jacobian_properties(assembly, ipoint, prescribed_conditions, point_masses)

Calculate/extract the point properties needed to calculate the mass matrix jacobian entries corresponding to a point for a constant mass matrix system

source
GXBeam.insert_static_point_jacobians!Function
insert_static_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
    resultants)

Insert the jacobian entries corresponding to a point for a static analysis into the system jacobian matrix.

source
GXBeam.insert_steady_point_jacobians!Function
insert_steady_point_jacobians!(jacob, indices, force_scaling, ipoint,
    properties, resultants, velocities)

Insert the jacobian entries corresponding to a point for a steady state analysis into the system jacobian matrix.

source
GXBeam.insert_initial_point_jacobians!Function
insert_initial_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
    resultants, velocities)

Insert the jacobian entries corresponding to a point for the initialization of a time domain analysis into the system jacobian matrix.

source
GXBeam.insert_dynamic_point_jacobians!Function
insert_dynamic_point_jacobians!(jacob, indices, force_scaling, ipoint,
    properties, resultants, velocities)

Insert the jacobian entries corresponding to a point for a steady state analysis into the system jacobian matrix.

source
GXBeam.insert_expanded_steady_point_jacobians!Function
insert_expanded_steady_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
    resultants, velocities)

Insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix for a constant mass matrix system.

source
GXBeam.insert_expanded_dynamic_point_jacobians!Function
insert_expanded_dynamic_point_jacobians!(jacob, indices, force_scaling, ipoint, properties,
    resultants, velocities)

Insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix for a constant mass matrix system.

source
GXBeam.insert_mass_matrix_point_jacobians!Function
insert_mass_matrix_point_jacobians!(jacob, gamma, indices, two_dimensional, force_scaling, ipoint,
    properties, resultants, velocities)

Insert the mass matrix jacobian entries corresponding to a point into the system jacobian matrix.

source
GXBeam.insert_expanded_mass_matrix_point_jacobians!Function
insert_expanded_mass_matrix_point_jacobians!(jacob, gamma, indices, two_dimensional, force_scaling, ipoint,
    properties, resultants, velocities)

Insert the mass matrix jacobian entries corresponding to a point into the system jacobian matrix.

source
GXBeam.static_point_jacobian!Function
static_point_jacobian!(jacob, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity)

Calculate and insert the jacobian entries corresponding to a point for a static analysis into the system jacobian matrix.

source
GXBeam.steady_point_jacobian!Function
steady_point_jacobian!(jacob, x, indices, force_scaling, assembly,
    ipoint, prescribed_conditions, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Calculate and insert the jacobian entries corresponding to a point for a steady state analysis into the system jacobian matrix.

source
GXBeam.initial_point_jacobian!Function
initial_point_jacobian!(jacob, x, indices, rate_vars,
    force_scaling, assembly, ipoint, prescribed_conditions, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration,
    u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate and insert the jacobian entries corresponding to a point for the initialization of a time domain analysis into the system jacobian matrix.

source
GXBeam.newmark_point_jacobian!Function
newmark_point_jacobian!(jacob, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity,
    udot_init, θdot_init, Vdot_init, Ωdot_init, dt))

Calculate and insert the jacobian entries corresponding to a point for a Newmark scheme time marching analysis into the system jacobian matrix.

source
GXBeam.dynamic_point_jacobian!Function
dynamic_point_jacobian!(jacob, dx, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)

Calculate and insert the jacobian entries corresponding to a point for a dynamic analysis into the system jacobian matrix.

source
GXBeam.expanded_steady_point_jacobian!Function
expanded_steady_point_jacobian!(jacob, x, indices, force_scaling,
    assembly, ipoint, prescribed_conditions, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Calculate and insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix.

source
GXBeam.expanded_dynamic_point_jacobian!Function
expanded_dynamic_point_jacobian!(jacob, dx, x, indices, force_scaling, assembly, ipoint,
    prescribed_conditions, point_masses, gravity, linear_velocity, angular_velocity)

Calculate and insert the jacobian entries corresponding to a point for a constant mass matrix system into the system jacobian matrix.

source
GXBeam.mass_matrix_point_jacobian!Function
mass_matrix_point_jacobian!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
    ipoint, prescribed_conditions)

Calculate and insert the mass_matrix jacobian entries corresponding to a point into the system jacobian matrix.

source
GXBeam.expanded_mass_matrix_point_jacobian!Function
expanded_mass_matrix_point_jacobian!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
    ipoint, prescribed_conditions, point_masses)

Calculate and insert the mass_matrix jacobian entries corresponding to a point into the system jacobian matrix for a constant mass matrix system

source

Elements

GXBeam.ElementType
Element{TF}

Composite type that defines a beam element's properties

Fields

  • L: Beam element length
  • x: Beam element location
  • compliance: Beam element compliance matrix
  • mass: Beam element mass matrix
  • Cab: Transformation matrix from the undeformed beam element frame to the body frame
  • mu: Beam element damping coefficients
source
GXBeam.element_loadsFunction
element_loads(x, ielem, icol_elem, force_scaling)

Extract the internal loads (F, M) for a beam element from the state variable vector. These loads are expressed in the deformed element frame.

source
GXBeam.expanded_element_loadsFunction
expanded_element_loads(x, ielem, icol_elem, force_scaling)

Extract the internal loads of a beam element at its endpoints (F1, M1, F2, M2) from the state variable vector for a constant mass matrix system.

source
GXBeam.expanded_element_velocitiesFunction
expanded_element_velocities(x, ielem, icol_elem)

Extract the velocities of a beam element from the state variable vector for a constant mass matrix system.

source
GXBeam.static_element_propertiesFunction
static_element_properties(x, indices, force_scaling, assembly, ielem,
    prescribed_conditions, gravity)

Calculate/extract the element properties needed to construct the residual for a static analysis

source
GXBeam.steady_element_propertiesFunction
steady_element_properties(x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
    linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))

Calculate/extract the element properties needed to construct the residual for a steady state analysis

source
GXBeam.initial_element_propertiesFunction
initial_element_properties(x, indices, rate_vars, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate/extract the element properties needed to construct the residual for a time-domain analysis initialization

source
GXBeam.newmark_element_propertiesFunction
newmark_element_properties(x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
    Vdot_init, Ωdot_init, dt)

Calculate/extract the element properties needed to construct the residual for a newmark- scheme time stepping analysis

source
GXBeam.dynamic_element_propertiesFunction
dynamic_element_properties(dx, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity)

Calculate/extract the element properties needed to construct the residual for a dynamic analysis

source
GXBeam.expanded_steady_element_propertiesFunction
expanded_steady_element_properties(x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity,
    linear_acceleration=(@SVector zeros(3)), angular_acceleration=(@SVector zeros(3)))

Calculate/extract the element properties needed to construct the residual for a constant mass matrix system

source
GXBeam.expanded_dynamic_element_propertiesFunction
expanded_dynamic_element_properties(dx, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, gravity, linear_velocity, angular_velocity)

Calculate/extract the element properties needed to construct the residual for a constant mass matrix system

source
GXBeam.static_element_resultantsFunction
static_element_resultants(properties, distributed_loads, ielem)

Calculate the resultant loads applied at each end of a beam element for a static analysis.

source
GXBeam.dynamic_element_resultantsFunction
dynamic_element_resultants(properties, distributed_loads, ielem)

Calculate the resultant loads applied at each end of a beam element for a steady state analysis.

source
GXBeam.insert_element_residuals!Function
insert_element_residuals!(resid, indices, force_scaling, assembly, ielem,
    compatibility, resultants)

Insert the residual entries corresponding to a beam element into the system residual vector.

source
GXBeam.insert_expanded_element_residuals!Function
insert_expanded_element_residuals!(resid, indices, force_scaling, assembly, ielem,
    compatibility, velocities, equilibrium, resultants)

Insert the residual entries corresponding to a beam element into the system residual vector for a constant mass matrix system.

source
GXBeam.static_element_residual!Function
static_element_residual!(resid, x, indices, force_scaling, assembly, ielem,
    prescribed_conditions, distributed_loads, gravity)

Calculate and insert the residual entries corresponding to a beam element for a static analysis into the system residual vector.

source
GXBeam.steady_element_residual!Function
steady_element_residual!(resid, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, distributed_loads, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Calculate and insert the residual entries corresponding to a beam element for a steady state analysis into the system residual vector.

source
GXBeam.initial_element_residual!Function
initial_element_residual!(resid, x, indices, rate_vars,
    force_scaling, structural_damping, assembly, ielem, prescribed_conditions,
    distributed_loads, gravity, linear_velocity, angular_velocity,
    linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate and insert the residual entries corresponding to a beam element for the initialization of a time domain simulation into the system residual vector.

source
GXBeam.newmark_element_residual!Function
newmark_element_residual!(resid, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, distributed_loads, gravity,
    linear_velocity, angular_velocity, Vdot_init, Ωdot_init, dt)

Calculate and insert the residual entries corresponding to a beam element for a newmark-scheme time marching analysis into the system residual vector.

source
GXBeam.dynamic_element_residual!Function
dynamic_element_residual!(resid, dx, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, distributed_loads, gravity,
    linear_velocity, angular_velocity)

Calculate and insert the residual entries corresponding to a beam element for a dynamic analysis into the system residual vector.

source
GXBeam.expanded_steady_element_residual!Function
expanded_steady_element_residual!(resid, x, indices, force_scaling, structural_damping,
    assembly, ielem, prescribed_conditions, distributed_loads, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Calculate and insert the residual entries corresponding to a beam element for a constant mass matrix system into the system residual vector.

source
GXBeam.expanded_dynamic_element_residual!Function
expanded_dynamic_element_residual!(resid, dx, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity)

Calculate and insert the residual entries corresponding to a beam element for a constant mass matrix system into the system residual vector.

source
GXBeam.static_element_jacobian_propertiesFunction
static_element_jacobian_properties(properties, x, indices, force_scaling,
    assembly, ielem, prescribed_conditions, gravity)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a beam element for a static analysis

source
GXBeam.steady_element_jacobian_propertiesFunction
steady_element_jacobian_properties(properties, x, indices,
    force_scaling, structural_damping, assembly, ielem, prescribed_conditions,
    gravity)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a steady state analysis

source
GXBeam.initial_element_jacobian_propertiesFunction
initial_element_jacobian_properties(properties, x, indices, rate_vars,
    force_scaling, structural_damping, assembly, ielem, prescribed_conditions, gravity,
    u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a Newmark scheme time marching analysis

source
GXBeam.newmark_element_jacobian_propertiesFunction
newmark_element_jacobian_properties(properties, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, gravity,
    Vdot_init, Ωdot_init, dt)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a Newmark scheme time marching analysis

source
GXBeam.dynamic_element_jacobian_propertiesFunction
dynamic_element_jacobian_properties(properties, dx, x, indices,
    force_scaling, assembly, ielem, prescribed_conditions, gravity)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a dynamic analysis

source
GXBeam.expanded_steady_element_jacobian_propertiesFunction
expanded_element_jacobian_properties(properties, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, gravity)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a steady state analysis

source
GXBeam.expanded_dynamic_element_jacobian_propertiesFunction
expanded_dynamic_element_jacobian_properties(properties, dx, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, gravity)

Calculate/extract the element properties needed to calculate the jacobian entries corresponding to a element for a dynamic analysis

source
GXBeam.mass_matrix_element_jacobian_propertiesFunction
mass_matrix_element_jacobian_properties(x, indices, force_scaling,
assembly, ielem, prescribed_conditions)

Calculate/extract the element properties needed to calculate the mass matrix jacobian entries corresponding to a beam element

source
GXBeam.expanded_mass_matrix_element_jacobian_propertiesFunction
expanded_mass_matrix_element_jacobian_properties(assembly, ielem, prescribed_conditions)

Calculate/extract the element properties needed to calculate the mass matrix jacobian entries corresponding to a beam element for a constant mass matrix system

source
GXBeam.static_element_resultant_jacobiansFunction
static_element_resultant_jacobians(properties, distributed_loads, ielem)

Calculate the jacobians for the resultant loads applied at each end of a beam element for a static analysis.

source
GXBeam.steady_element_resultant_jacobiansFunction
steady_element_resultant_jacobians(properties, distributed_loads, ielem)

Calculate the jacobians for the resultant loads applied at each end of a beam element for a steady state analysis.

source
GXBeam.initial_element_resultant_jacobiansFunction
initial_element_resultant_jacobians(properties, distributed_loads, ielem)

Calculate the jacobians for the resultant loads applied at each end of V beam element for the initialization of a time domain analysis.

source
GXBeam.dynamic_element_resultant_jacobiansFunction
dynamic_element_resultant_jacobians(properties, distributed_loads, ielem)

Calculate the jacobians for the resultant loads applied at each end of a beam element for a dynamic analysis.

source
GXBeam.static_element_jacobian!Function
static_element_jacobian!(jacob, x, indices, force_scaling,
    assembly, ielem, prescribed_conditions, distributed_loads, gravity)

Calculate and insert the jacobian entries corresponding to a beam element for a static analysis into the system jacobian matrix.

source
GXBeam.steady_element_jacobian!Function
steady_element_jacobian!(jacob, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Calculate and insert the jacobian entries corresponding to a beam element for a steady state analysis into the system jacobian matrix.

source
GXBeam.initial_element_jacobian!Function
initial_element_jacobian!(jacob, x, indices, rate_vars, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, ub, θb, vb, ωb, ab, αb, u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Calculate and insert the jacobian entries corresponding to a beam element for the initialization of a time domain analysis into the system jacobian matrix.

source
GXBeam.newmark_element_jacobian!Function
newmark_element_jacobian!(jacob, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity, Vdot_init, Ωdot_init, dt)

Calculate and insert the jacobian entries corresponding to a beam element for a Newmark-scheme time marching analysis into the system jacobian matrix.

source
GXBeam.dynamic_element_jacobian!Function
dynamic_element_jacobian!(jacob, dx, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity)

Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.

source
GXBeam.expanded_steady_element_jacobian!Function
expanded_steady_element_jacobian!(jacob, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity, linear_acceleration,
    angular_acceleration)

Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.

source
GXBeam.expanded_dynamic_element_jacobian!Function
expanded_dynamic_element_jacobian!(jacob, dx, x, indices, force_scaling,
    structural_damping, assembly, ielem, prescribed_conditions, distributed_loads,
    gravity, linear_velocity, angular_velocity)

Calculate and insert the jacobian entries corresponding to a beam element for a dynamic analysis into the system jacobian matrix.

source
GXBeam.mass_matrix_element_jacobian!Function
mass_matrix_element_jacobian!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
    ielem, prescribed_conditions)

Calculate and insert the mass_matrix jacobian entries corresponding to a beam element into the system jacobian matrix.

source
GXBeam.expanded_mass_matrix_element_jacobian!Function
expanded_mass_matrix_element_jacobian!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
    ielem, prescribed_conditions)

Calculate and insert the mass_matrix jacobian entries corresponding to a beam element into the system jacobian matrix for a constant mass matrix system

source

System

GXBeam.SystemIndicesType
SystemIndices

Structure for holding indices for accessing the state variables and equations associated with each point and beam element in a system.

source
GXBeam.default_force_scalingFunction
default_force_scaling(assembly)

Defines a suitable default force scaling factor based on the nonzero elements of the compliance matrices in assembly.

source
GXBeam.curve_triadFunction
curve_triad(Cab, k, s)
curve_triad(Cab, kkt, ktilde, kn, s)

Return the transformation matrix at s along the length of the beam given the curvature vector k and the initial transformation matrix Cab.

source
GXBeam.curve_coordinatesFunction
curve_coordiantes(r, Cab, k, s)
curve_coordinates(r, Cab, kkt, ktilde, kn, s)

Return the coordinates at s along the length of the beam given the starting point r, initial transformation matrix Cab, and curvature vector k.

source
GXBeam.set_initial_state!Function
set_initial_state!([x=system.x,] system, rate_vars, assembly, state; kwargs...)

Set the state variables in x for an initial condition analysis to the values in state

source
GXBeam.static_system_residual!Function
static_system_residual!(resid, x, indices, two_dimensional, force_scaling,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity)

Populate the system residual vector resid for a static analysis

source
GXBeam.initial_system_residual!Function
initial_system_residual!(resid, x, indices, rate_vars1, rate_vars2,
    two_dimensional, force_scaling, structural_damping, assembly, prescribed_conditions,
    distributed_loads, point_masses, gravity, linear_velocity, angular_velocity,
    linear_acceleration, angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Populate the system residual vector resid for the initialization of a time domain simulation.

source
GXBeam.steady_system_residual!Function
steady_system_residual!(resid, x, indices, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, vb, ωb, ab, αb)

Populate the system residual vector resid for a steady state analysis

source
GXBeam.newmark_system_residual!Function
newmark_system_residual!(resid, x, indices, two_dimensional, force_scaling, structural_damping,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
    linear_velocity, angular_velocity, udot_init, θdot_init, Vdot_init, Ωdot_init, dt)

Populate the system residual vector resid for a Newmark scheme time marching analysis.

source
GXBeam.dynamic_system_residual!Function
dynamic_system_residual!(resid, dx, x, indices, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, linear_velocity, angular_velocity)

Populate the system residual vector resid for a general dynamic analysis.

source
GXBeam.expanded_steady_system_residual!Function
expanded_steady_system_residual!(resid, x, indices, two_dimensional, force_scaling, structural_damping,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
    linear_velocity, angular_velocity, linear_acceleration, angular_acceleration)

Populate the system residual vector resid for a constant mass matrix system.

source
GXBeam.expanded_dynamic_system_residual!Function
expanded_dynamic_system_residual!(resid, dx, x, indices, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, linear_velocity, angular_velocity)

Populate the system residual vector resid for a constant mass matrix system.

source
GXBeam.static_system_jacobian!Function
static_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity)

Populate the system jacobian matrix jacob for a static analysis

source
GXBeam.steady_system_jacobian!Function
steady_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, linear_velocity, angular_velocity, linear_acceleration,
    angular_acceleration)

Populate the system jacobian matrix jacob for a steady-state analysis

source
GXBeam.initial_system_jacobian!Function
initial_system_jacobian!(jacob, x, indices, rate_vars1, rate_vars2, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, linear_velocity, angular_velocity, linear_acceleration,
    angular_acceleration, u0, θ0, V0, Ω0, Vdot0, Ωdot0)

Populate the system jacobian matrix jacob for the initialization of a time domain simulation.

source
GXBeam.newmark_system_jacobian!Function
newmark_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling, structural_damping,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
    linear_velocity, angular_velocity, udot_init, θdot_init, Vdot_init, Ωdot_init, dt)

Populate the system jacobian matrix jacob for a Newmark scheme time marching analysis.

source
GXBeam.dynamic_system_jacobian!Function
dynamic_system_jacobian!(jacob, dx, x, indices, two_dimensional, force_scaling,
    structural_damping, assembly, prescribed_conditions, distributed_loads,
    point_masses, gravity, linear_velocity, angular_velocity)

Populate the system jacobian matrix jacob for a general dynamic analysis.

source
GXBeam.expanded_steady_system_jacobian!Function
expanded_steady_system_jacobian!(jacob, x, indices, two_dimensional, force_scaling, structural_damping,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
    ub_p, θb_p, vb_p, ωb_p, ab_p, αb_p)

Populate the system jacobian matrix jacob for a general dynamic analysis with a constant mass matrix system.

source
GXBeam.expanded_dynamic_system_jacobian!Function
expanded_dynamic_system_jacobian!(jacob, dx, x, indices, two_dimensional, force_scaling, structural_damping,
    assembly, prescribed_conditions, distributed_loads, point_masses, gravity,
    linear_velocity, angular_velocity)

Populate the system jacobian matrix jacob for a general dynamic analysis with a constant mass matrix system.

source
GXBeam.system_mass_matrix!Function
system_mass_matrix!(jacob, x, indices, two_dimensional, force_scaling,  assembly,
    prescribed_conditions, point_masses)

Calculate the jacobian of the residual expressions with respect to the state rates.

source
system_mass_matrix!(jacob, gamma, x, indices, two_dimensional, force_scaling, assembly,
    prescribed_conditions, point_masses)

Calculate the jacobian of the residual expressions with respect to the state rates and add the result multiplied by gamma to jacob.

source
GXBeam.expanded_system_mass_matrixFunction
expanded_system_mass_matrix(system, assembly;
    two_dimensional = false,
    prescribed_conditions=Dict{Int, PrescribedConditions}(),
    point_masses=Dict{Int, PointMass}())

Calculate the jacobian of the residual expressions with respect to the state rates for a constant mass matrix system.

source
GXBeam.expanded_system_mass_matrix!Function
expanded_system_mass_matrix!(jacob, indices, two_dimensional, force_scaling,  assembly, prescribed_conditions,
    point_masses)

Calculate the jacobian of the residual expressions with respect to the state rates.

source
expanded_system_mass_matrix!(jacob, gamma, indices, two_dimensional, force_scaling, assembly,
    prescribed_conditions, point_masses)

Calculate the jacobian of the residual expressions with respect to the state rates and add the result multiplied by gamma to jacob.

source

Analysis

GXBeam.initialize_system!Function
initialize_system!(system, assembly, tvec; kwargs...)

Pre-allocate the history for a stepped time-domain analysis and conduct the initial condition analysis.

source
GXBeam.step_system!Function
step_system!(history, system, assembly, tvec; kwargs...)

A similar function to timedomainanalysis, but instead history is already allocated and it only takes one step with the system.

Inputs

  • system::DynamicSystem : The dynamic system to be stepped.
  • paug: The augmented parameter vector.
  • x: The current state vector.
  • constants: The constants tuple.
  • initial_state: The initial state of the system.
  • assembly: The assembly of the system.
  • tvec: The time vector.
  • i: The current time step index.
  • prescribed_conditions: A dictionary of prescribed conditions.
  • distributed_loads: A dictionary of distributed loads.
  • point_masses: A dictionary of point masses.
  • linear_velocity: The linear velocity of the body frame.
  • angular_velocity: The angular velocity of the body frame.
  • linear_acceleration: The linear acceleration of the body frame.
  • angular_acceleration: The angular acceleration of the body frame.
  • gravity: The gravity vector in the body frame.
  • structural_damping: A flag indicating whether to enable structural damping.
  • linear: A flag indicating whether a linear analysis should be performed.
  • show_trace: A flag indicating whether to display the solution progress.
  • update_linearization: A flag indicating whether to update the linearization state variables for a linear analysis with the instantaneous state variables.
  • xpfunc: A function that returns parameters as a function of the state variables.
  • pfunc: A function that returns parameters as a function of time.
  • p: Sensitivity parameters.
  • converged: A reference to a boolean indicating whether the solution converged.
source
GXBeam.take_stepFunction
take_step(x, dx, system, assembly, t, tprev, 
prescribed_conditions, distributed_loads, point_masses,
gravity, linear_velocity, angular_velocity,
xpfunc, pfunc, p,
structural_damping, two_dimensional)

Take a single implicit Euler time step for the given system.

source
GXBeam.simulateFunction
simulate(assembly, tvec; prescribed_conditions, distributed_loads, point_masses, gravity, linear_velocity, angular_velocity,
xpfunc, pfunc, p, structural_damping, two_dimensional, verbose)

Use the takestep function to run a timedomain_analysis.

WARNING: The implicit Euler method is not as stable as the Newmark-beta method. You may experience instability for stiff systems or large time steps.

source

Section

GXBeam.combine_halfsFunction

given nodes/elements for upper surface and lower surface separately, combine into one set making sure to reuse the common nodes that occur at the LE/TE

source
GXBeam.insertpointFunction

add new point xn into existing airfoil coordinates x, y. Used to ensure breakpoints are in airfoil mesh. assumes x is already sorted (and thus represents either upper or lower half)

source
GXBeam.tangentialFunction

determine tangential direction for each point on airfoil (as a normal vector with magnitudes tx and ty). keyword signifies whether this is the upper or lower surface

source
GXBeam.element_orientationFunction
element_orientation(nodes)

Calculate the orientation beta of a single element using its nodes. Return cos(beta) and sin(beta)

source
GXBeam.addwebsFunction

add the nodes and elements for the webs. given the x locations (by idx) where the webs start and the number of grid points in the webs.

source
GXBeam.web_intersectionsFunction

determine where web intersects the inner surface curves. creates a new discretization where points are added to line up with the vertical mesh in the web. also returns the x vector index where the web mesh should start

source
GXBeam.parseairfoilFunction

airfoil coordinates xaf, yaf are assumed to start at trailing edge and traverse counterclockwise back to trailing edge. (thus upper surface first) coordinates are also assumed to be chord-normalized so that xaf is between 0 and 1. xbreak is a vector defining locations that need to be inserted into airfoil coordinates (for upper and lower surfaces) ds is a desired spacing for resampling (set ds=nothing to not resample) Alternatively can use vector nseg to define number of points in each interval between xbreak (see resample())

Method separates airfoil into upper and lower surfaces, adds in break points, and optionally resamples.

source
GXBeam.node2idxFunction
node2idx(nodenums)

Convenience function to map node numbers to locations in global matrix.

source
GXBeam.reorderFunction
reorder(K)

Reorder stiffness or compliance matrix from internal order to GXBeam order

source
GXBeam.elementQFunction
elementQ(material, theta, cbeta, sbeta)

Calculate element constitutive matrix Q accounting for fiber orientation theta and element orientation beta (where cbeta = cos(beta) and sbeta=sin(beta))

source
GXBeam.resampleFunction

resample airfoil coordinates x, y with linear interpolation, making sure to exactly hit the points in vector xbreak (e.g., [0, 0.2, 0.4, 0.7, 1.0]). ds defines the target size between coordinates (rounded to nearnest number of grid points). assumes x is already sorted (and thus represents either upper or lower half or airfoil coordinates) points in xbreak do not have to exist in the original x vector alternatively you can define nseg, a vector of size: length(xbreak) - 1, defining how many points to use in each interval (e.g., [15, 20, 40, 30]) If nseg is provided then ds is ignored. During optimization either keep mesh fixed, or if changes are needed, then nseg should be used rather than ds to avoid discrete changes.

source
GXBeam.stiffnessFunction
stiffness(material)

Construct constitutive matrix Q for the specified material (with the internal ordering).

source
GXBeam.preprocess_layersFunction

convert segments to all have the same number of layers for ease in meshing. Overall definition remains consistent, just break up some thicker layers into multiple thinner layers (with same material and orientation properties)

source
GXBeam.linearsolveFunction
linearsolve(A, B; kwargs...)

Linear solve which is overloaded with analytic derivatives. Returns the factorized matrix Afact and the result of the linear solve X.

source
GXBeam.element_integrandFunction
element_integrand(ksi, eta, element, nodes)

Compute the integrand for a single element with a given ksi, eta.

source
GXBeam.te_inner_intersectionFunction

given curves for inner surfaces, find where they intersect (if they intersect), for handling trailing-edge mesh returns index of upper and lower surface where the trailing-edge-mesh should start (before intersection). also returns the x, y location of intersection.

source
GXBeam.mesh_areaFunction
mesh_area(nodes, elements)

Compute total area of the mesh.

Inputs

  • nodes::Vector{Node{TF}}: all the nodes in the mesh
  • elements::Vector{MeshElement{TF}}: all the elements in the mesh

Returns

  • A::TF: total area of the mesh
source
GXBeam.mesh_cylinderFunction
mesh_cylinder(R, thickness, material)

A simple function to mesh a circular cross-section.

Inputs

  • r::Vector{Float} - A vector of all the radial location of nodes. Assumed from least to greatest.
  • materials::Vector{Material} - All of the materials used in the cross section.
  • material_idx::Vector{Int} - The index of the material in the materials vector for each layer.

Assumed in order of outermost layer to inner most layer. Note: This is opposite to the order of r.

  • plyangles::Vector{Float} - The ply angle of a given radial element (same order as material_idx).
  • nt::Int - Number of tangential points (How many elements around the circle).
source

Index