Library

Public API

Creating an Assembly

GXBeam.curve_lengthFunction
curve_length(start, stop, curvature)

Calculate the length of a curve given its endpoints and its curvature vector

source
GXBeam.discretize_beamFunction
discretize_beam(L, start, discretization; frame = Matrix(I,3,3)), curvature = zeros(3))

Discretize a beam according to the discretization provided in discretization given the beam length (L), and starting point (start).

Return the lengths, endpoints, midpoints, and rotation matrices of the beam elements.

Arguments

  • L: Beam length
  • start: Beam starting point
  • discretization: May be either an integer, representing the number of elements that the beam should be discretized into, or a vector containing the normalized endpoints of each beam element, where 0 is the beginning of the beam and 1 is the end of the beam.
  • frame: 3x3 beam rotation matrix which transforms from the local beam coordinate frame at the start of the beam to the global coordinate frame.
  • curvature: curvature vector
source
GXBeam.AssemblyMethod
Assembly(points, start, stop; kwargs...)

Construct an assembly of connected nonlinear beam elements for analysis. Beam lengths and midpoints may be manually specified in case beam elements are curved rather than straight.

Arguments

  • points: Array of all beam element endpoints
  • start: Array containing point indices where each beam element starts
  • stop: Array containing point indices where each beam element stops

Keyword Arguments

  • stiffness: Array of (6 x 6) stiffness matrices for each beam element, alternative to providing compliance
  • compliance: Array of (6 x 6) compliance matrices for each beam element, defaults to zeros(6,6) for each beam element
  • mass: Array of (6 x 6) mass matrices for each beam element, alternative to providing minv
  • minv: Array of (6 x 6) mass matrices for each beam element, defaults to the identity matrix for each beam element
  • frames: Array of (3 x 3) rotation matrices for each beam element (to transform to the global frame), defaults to the identity matrix for each beam element
  • lengths: Array containing the length of each beam, defaults to the distance between beam endpoints
  • midpoints: Array containing the midpoint of each beam element, defaults to the average of the beam element endpoints
source

Defining Distributed Loads

GXBeam.DistributedLoadsMethod
DistributedLoads(assembly, ibeam[, dt]; kwargs...)

Integrates the specified distributed loads over the element for each time step.

Arguments

  • assembly: The beam element assembly
  • ibeam: The index of the beam element which the distributed load is assigned to
  • dt: Time step size. If omitted a single time step is assumed and specified functions become a function of s only.
  • s1 = 0.0: Start of beam element (used for integrating the distributed loads)
  • s2 = 1.0: End of beam element (used for integrating the distributed loads)
  • nstep: The total length of the time vector
  • method = (f, a, b) -> gauss_quadrature(f, a, b): Method which integrates function f from a to b. Defaults to the Gauss-Legendre quadrature with 4 points on each element.
  • fx = (s, t) -> 0.0: Distributed non-follower force on beam element in x-direction
  • fy = (s, t) -> 0.0: Distributed non-follower force on beam element in y-direction
  • fz = (s, t) -> 0.0: Distributed non-follower force on beam element in z-direction
  • mx = (s, t) -> 0.0: Distributed non-follower moment on beam element in x-direction
  • my = (s, t) -> 0.0: Distributed non-follower moment on beam element in y-direction
  • mz = (s, t) -> 0.0: Distributed non-follower moment on beam element in z-direction
  • fx_follower = (s, t) -> 0.0: Distributed follower force on beam element in x-direction
  • fy_follower = (s, t) -> 0.0: Distributed follower force on beam element in y-direction
  • fz_follower = (s, t) -> 0.0: Distributed follower force on beam element in z-direction
  • mx_follower = (s, t) -> 0.0: Distributed follower moment on beam element in x-direction
  • my_follower = (s, t) -> 0.0: Distributed follower moment on beam element in y-direction
  • mz_follower = (s, t) -> 0.0: Distributed follower moment on beam element in z-direction
source

Defining Prescribed Conditions

GXBeam.PrescribedConditionsMethod
PrescribedConditions(dt=0.0; kwargs...)

Construct an object of type PrescribedConditions which stores the prescribed conditions for a point at each time step.

Prescribed conditions may be assigned as either a scalar parameter or as a function of time.

Prescribed Wiener-Milenkovic parameters must satisfy the following inequality: sqrt(thetax^2 + thetay^2 + theta_z^2) <= 4. Note that this restriction still allows all possible rotations to be represented.

Arguments

  • dt: Time step size.
  • nstep: The total length of the time vector
  • ux: Prescribed x-direction displacement of the point
  • uy: Prescribed y-direction displacement of the point
  • uz: Prescribed z-direction displacement of the point
  • theta_x: Prescribed first Wiener-Milenkovic parameter of the point
  • theta_y: Prescribed second Wiener-Milenkovic parameter of the point
  • theta_z: Prescribed third Wiener-Milenkovic parameter of the point
  • Fx: Prescribed force in x-direction applied on the point
  • Fy: Prescribed force in y-direction applied on the point
  • Fz: Prescribed force in z-direction applied on the point
  • Mx: Prescribed moment about x-axis applied on the point
  • My: Prescribed moment about y-axis applied on the point
  • Mz: Prescribed moment about z-axis applied on the point
  • Fx_follower: Prescribed follower force in x-direction applied on the point
  • Fy_follower: Prescribed follower force in y-direction applied on the point
  • Fz_follower: Prescribed follower force in z-direction applied on the point
  • Mx_follower: Prescribed follower moment about x-axis applied on the point
  • My_follower: Prescribed follower moment about y-axis applied on the point
  • Mz_follower: Prescribed follower moment about z-axis applied on the point
source

Pre-Initializing Memory for an Analysis

GXBeam.SystemMethod
System([TF=eltype(assembly),] assembly, points, static)

Initialize an object of type System which stores the system state, residual vector, current time function values,and jacobian matrices as well as pointers to be able to access their contents.

Arguments:

  • TF: (optional) Used to specify floating point type used by resulting System object
  • assembly: Assembly of rigidly connected nonlinear beam elements
  • points: Point indices which should be preserved in the system of equations. All points with prescribed conditions should be included.
  • static: Flag indicating whether system matrices will be used for static simulations
source

Performing an Analysis

GXBeam.static_analysisFunction
static_analysis(assembly; kwargs...)

Perform a static analysis of the system of nonlinear beams contained in assembly. Return the resulting system and a flag indicating whether the iteration procedure converged.

Keyword Arguments

  • prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}(): Dictionary holding PrescribedConditions composite types for the points in keys(prescribed_conditions)
  • distributed_loads = Dict{Int,DistributedLoads{Float64}}(): Dictionary holding DistributedLoads composite types for the beam elements in keys(distributed_loads)
  • linear = false: Set to true for a linear analysis
  • method = :newton: Method (as defined in NLsolve) to solve nonlinear system of equations
  • linesearch = LineSearches.BackTracking(maxstep=1e6): Line search used to solve nonlinear system of equations
  • ftol = 1e-9: tolerance for solving nonlinear system of equations
  • iterations = 1000: maximum iterations for solving the nonlinear system of equations
  • nstep = 1: Number of time steps. May be used in conjunction with time varying prescribed conditions and distributed loads to gradually increase displacements/loads.
source
GXBeam.static_analysis!Function
static_analysis!(system, assembly; kwargs...)

Pre-allocated version of static_analysis. Uses the state variables stored in system as an initial guess for iterating.

source
GXBeam.steady_state_analysisFunction
steady_state_analysis(assembly; kwargs...)

Perform a steady-state analysis for the system of nonlinear beams contained in assembly. Return the resulting system and a flag indicating whether the iteration procedure converged.

Keyword Arguments

  • prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}(): Dictionary holding PrescribedConditions composite types for the points in keys(prescribed_conditions)
  • distributed_loads = Dict{Int,DistributedLoads{Float64}}(): Dictionary holding DistributedLoads composite types for the beam elements in keys(distributed_loads)
  • linear = false: Set to true for a linear analysis
  • method = :newton: Method (as defined in NLsolve) to solve nonlinear system of equations
  • linesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6): Line search used to solve nonlinear system of equations
  • ftol = 1e-9: tolerance for solving nonlinear system of equations
  • iterations = 1000: maximum iterations for solving the nonlinear system of equations
  • nstep = 1: Number of time steps. May be used in conjunction with time varying prescribed conditions, distributed loads, and global motion to gradually increase displacements/loads.
  • origin = zeros(3): Global frame origin
  • linear_velocity = fill(zeros(3), nstep): Global frame linear velocity for each time step.
  • angular_velocity = fill(zeros(3), nstep): Global frame angular velocity for each time step
source
GXBeam.steady_state_analysis!Function
steady_state_analysis!(system, assembly; kwargs...)

Pre-allocated version of steady_state_analysis. Uses the state variables stored in system as an initial guess for iterating.

source
GXBeam.eigenvalue_analysisFunction
eigenvalue_analysis(assembly; kwargs...)

Compute the eigenvalues and eigenvectors of the system of nonlinear beams contained in assembly. Return the modified system, eigenvalues, eigenvectors, and a convergence flag indicating whether the corresponding steady-state analysis converged.

Keyword Arguments

  • prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}(): Dictionary holding PrescribedConditions composite types for the points in keys(prescribed_conditions)
  • distributed_loads = Dict{Int,DistributedLoads{Float64}}(): Dictionary holding DistributedLoads composite types for the beam elements in keys(distributed_loads)
  • linear = false: Set to true for a linear analysis
  • method = :newton: Method (as defined in NLsolve) to solve nonlinear system of equations
  • linesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6): Line search used to solve nonlinear system of equations
  • ftol = 1e-9: tolerance for solving nonlinear system of equations
  • iterations = 1000: maximum iterations for solving the nonlinear system of equations
  • nstep = 1: Number of time steps. May be used in conjunction with time varying prescribed conditions, distributed loads, and global motion to gradually increase displacements/loads.
  • origin = zeros(3): Global frame origin
  • linear_velocity = fill(zeros(3), nstep): Global frame linear velocity for each time step.
  • angular_velocity = fill(zeros(3), nstep): Global frame angular velocity for each time step
  • nev = 6: Number of eigenvalues to compute
source
GXBeam.eigenvalue_analysis!Function
eigenvalue_analysis!(system, assembly; kwargs...)

Pre-allocated version of eigenvalue_analysis. Uses the state variables stored in system as an initial guess for iterating to find the steady state solution.

source
GXBeam.time_domain_analysisFunction
time_domain_analysis(assembly, dt; kwargs...)

Perform a time-domain analysis for the system of nonlinear beams contained in assembly. Return the final system, a post-processed solution history, and a convergence flag indicating whether the iterations converged for each time step.

Keyword Arguments

  • prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}(): Dictionary holding PrescribedConditions composite types for the points in keys(prescribed_conditions)
  • distributed_loads = Dict{Int,DistributedLoads{Float64}}(): Dictionary holding DistributedLoads composite types for the beam elements in keys(distributed_loads)
  • linear = false: Set to true for a linear analysis
  • method = :newton: Method (as defined in NLsolve) to solve nonlinear system of equations
  • linesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6): Line search used to solve nonlinear system of equations
  • ftol = 1e-9: tolerance for solving nonlinear system of equations
  • iterations = 1000: maximum iterations for solving the nonlinear system of equations
  • nstep = 1: The total length of the time vector
  • origin = zeros(3): Global frame origin
  • linear_velocity = fill(zeros(3), nstep): Global frame linear velocity for each time step.
  • angular_velocity = fill(zeros(3), nstep): Global frame angular velocity for each time step.
  • u0=fill(zeros(3), length(assembly.elements)): initial displacment of each beam element,
  • theta0=fill(zeros(3), length(assembly.elements)): initial angular displacement of each beam element,
  • udot0=fill(zeros(3), length(assembly.elements)): initial time derivative with respect to u
  • thetadot0=fill(zeros(3), length(assembly.elements)): initial time derivative with respect to theta
  • save=1:nstep: Steps at which to save the time history
source
GXBeam.time_domain_analysis!Function
time_domain_analysis!(system, assembly, dt; kwargs...)

Pre-allocated version of time_domain_analysis. Uses the state variables stored in system as an initial guess for iterating.

source

Post-Processing

GXBeam.AssemblyStateMethod
AssemblyState(system, assembly, x = system.x;
    prescribed_conditions = Dict{Int,PrescribedConditions{Float64}}())

Post-process the system state given the solution vector x. Return an object of type AssemblyState that defines the state of the assembly for the time step.

If prescribed_conditions is not provided, all point state variables are assumed to be displacements/rotations, rather than their actual identities as used in the analysis.

source
GXBeam.AssemblyStateType
AssemblyState{TF, TP<:AbstractVector{PointState{TF}}, TE<:AbstractVector{ElementState{TF}}}

Struct for storing state variables for the points and elements in an assembly.

Fields:

  • points::TP: Array of PointStates for each point in the assembly
  • elements::TE: Array of ElementStates for each element in the assembly
source
GXBeam.left_eigenvectorsFunction
left_eigenvectors(system, λ, V)
left_eigenvectors(K, M, λ, V)

Compute the left eigenvector matrix U for the system using inverse power iteration given the eigenvalues λ and the corresponding right eigenvector matrix V.

The complex conjugate of each left eigenvector is stored in each row of the matrix U

Left and right eigenvectors satisfy the following M-orthogonality condition:

  • u'Mv = 1 if u and v correspond to the same eigenvalue
  • u'Mv = 0 if u and v correspond to different eigenvalues

This means that UMV = I

This function assumes that system has not been modified since the eigenvalues and right eigenvectors were computed.

source
GXBeam.correlate_eigenmodesFunction
correlate_eigenmodes(C)

Return the permutation and the associated corruption index vector which associates eigenmodes from the current iteration with those of the previous iteration given the correlation matrix C.

The correlation matrix can take one of the following forms (in order of preference):

  • C = U_p*M*V
  • C = U*M_p*V_p
  • C = V_p'*V
  • C = V'*V_p

where U is a matrix of conjugated left eigenvectors, M is the system mass matrix, V is a matrix of right eigenvectors, and ()_p indicates a variable from the previous iteration.

Note that the following two forms of the correlation matrix seem to be significantly inferior to their counterparts listed above: C = U*M*V_p and C = U_p*M_p*V. This is likely due to the way in which the left eigenvector matrix is calculated.

The corruption index is the largest magnitude in a given row of C that was not chosen divided by the magnitude of the chosen eigenmode. It is most meaningful when using one of the forms of the correlation matrix that uses left eigenvectors since correct eigenmodes will have magnitudes close to 1 and incorrect eigenmodes will have magnitudes close to 0.

If the new mode number is already assigned, the next highest unassigned mode number is used. In this case a corruption index higher than 1 will be returned, otherwise the values of the corruption index will always be bounded by 0 and 1.

See "New Mode Tracking Methods in Aeroelastic Analysis" by Eldred, Vankayya, and Anderson.

source
GXBeam.wiener_milenkovicFunction
wiener_milenkovic(c)

Construct a Wiener-Milenkovic rotation matrix, given the three Wiener-Milenkovic parameters in c.

source
GXBeam.write_vtkFunction
write_vtk(name, assembly::Assembly; kwargs...)
write_vtk(name, assembly::Assembly, state::AssemblyState; kwargs...)
write_vtk(name, assembly::Assembly, history::Vector{<:AssemblyState}], dt; kwargs...)

Write the deformed geometry (and associated data) to a VTK file for visualization using ParaView.

The state parameter may be omitted to write the original geometry to a VTK file without any associated data.

If the solution time history is provided, the time step must also be provided

Keyword Arguments

  • scaling=1.0: Parameter to scale the deflections (only valid if state is provided)
  • metadata=Dict(): Dictionary of metadata for the file(s)
source
write_vtk(name, assembly::Assembly, [state::AssemblyState, ]λ::Number, eigenstate::AssemblyState;
scaling=1.0, mode_scaling=1.0, cycles=1, steps=100)

Write a series of files corresponding to the elastic motion of the assembly about the deformed state encoded in state defined by the eigenvalue λ and the eigenvector encoded in eigenstate over the time period specified by time.

The steady-state deflections can be scaled with scaling and the eigenmode deflections can be scaled using mode_scaling.

The current time is encoded in the metadata tag "time"

source

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 rotation matrix C given the three angular parameters in θ.

source
GXBeam.get_C_tFunction
get_C_t([C, ] θ, θ_t)

Calculate the derivative of the Wiener-Milenkovic rotation matrix C with respect to the scalar parameter t. θ_t is the derivative of the angular parameter θ with respect to t.

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

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

source
GXBeam.get_C_θdotFunction
get_C_θdot([C, ] θ)

Calculate the derivative of the time derivative of the Wiener-Milenkovic rotation matrix C with respect to each of the time derivatives of θ. Used for constructing the "mass" matrix for eigenvalue computations.

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 c.

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.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

Points

GXBeam.point_variablesFunction
point_variables(x, icol)
point_variables(x, icol, prescribed_conditions, istep)

Extract u, θ, F, M for the point described by the point state variables at icol in x after incorporating the prescribed conditions in prescribed_conditions

Note that the degrees of freedom that are not specified in prescribed_conditions are used as state variables (e.g. prescribing F[2] would mean u[2] = x[icol+1])

source
GXBeam.insert_point_residual!Function
insert_point_residual!(resid, irow_b, irow_p, u, θ, F, M, side)

Modify the equilibrium and constitutive equations to account for the point variables given by u, θ, F, M

If irowb != irowp, assume that the equilibrium equations have already been modified

Arguments

  • resid: Residual vector
  • irow_b: Row index of the first equilibrium/compatability equation for one side of the beam element
  • irow_p: Row index of the first equilibrium equation for the point
  • u: Displacement of the point
  • θ: Rotation of the point
  • F: External forces imposed on the point
  • M: External moments imposed on the point
  • side: Side of beam (-1 (left) or 1 (right))
source
GXBeam.point_residual!Function
point_residual!(resid, x, ipoint, assembly, prescribed_conditions, icol,
    irow_p, irow_beam1, irow_beam2)

Adds a points contributions to the residual vector

Arguments

  • resid: residual vector
  • x: current state vector
  • ipoint: index of point
  • assembly: assembly of interconnected beam elements
  • prescribed_conditions: dictionary of prescribed conditions
  • icol: starting index for the point's state variables
  • irow_p: Row index of the first equilibrium equation for the point
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
source
GXBeam.point_follower_jacobiansFunction
point_follower_jacobians(x, icol, prescribed_conditions)

Calculate the jacobians of the follower forces/moments with respect to θ

Arguments

  • x: Current state variable vector
  • icol: Row/Column index of the first state variable for the point
  • prescribed_conditions: Prescribed conditions for the point
source
GXBeam.insert_point_jacobian!Function
insert_point_jacobian!(jacob, irow_b, irow_p, icol, prescribed_conditions, side, F_θ, M_θ)

Modify the jacobian entries for the equilibrium and constitutive equations to account for the point variables at icol

If irowb != irowp, assume that the equilibrium equations have already been modified

Arguments

  • jacob: Jacobian of residual vector with respect to state vectors
  • irow_b: Row index of the first equilibrium/compatability equation for one side of the beam element
  • irow_p: Row index of the first equilibrium equation for the point
  • icol: Row/Column index of the first state variable for the point
  • prescribed_conditions: Prescribed force/displacement and moment/rotation on point
  • side: Side of beam (-1 (left) or 1 (right))
source
GXBeam.point_jacobian!Function
point_jacobian!(jacob, x, ipoint, assembly, prescribed_conditions, icol,
    irow_p, irow_beam1, irow_beam2)

Adds a points contributions to the residual vector

Arguments

  • jacob: residual vector
  • x: current state vector
  • ipoint: index of point
  • assembly: assembly of interconnected beam elements
  • prescribed_conditions: dictionary of prescribed conditions
  • istep: current time step
  • icol: starting index for the point's state variables
  • irow_p: Row index of the first equilibrium equation for the point
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
source
GXBeam.PointStateType
PointState

Holds the state variables for a point

Fields:

  • u: Displacement variables for the point (in the global coordinate frame)
  • theta: Wiener-Milenkovic rotational displacement variables for the point
  • F: Externally applied forces on the point (in the global coordinate frame)
  • M: Externally applied moments on the point (in the global coordinate frame)
source

Elements

GXBeam.ElementType
Element{TF}

Composite type that defines a beam element's properties

Fields

  • L: Length of the beam element
  • x: Location of the beam element (the midpoint of the beam element)
  • C11: Upper left portion of the beam element's compliance matrix
  • C12: Upper right portion of the beam element's compliance matrix
  • C22: Lower right portion of the beam element's compliance matrix
  • minv11: Upper left portion of the inverse of the beam element's mass matrix
  • minv12: Upper right portion of the inverse of the beam element's mass matrix
  • minv22: Lower right portion of the inverse of the beam element's mass matrix
  • Cab: Rotation matrix from the global frame to beam element's frame
source
GXBeam.element_strainFunction
element_strain(element, F, M)

Calculate the strain of a beam element given the resultant force and moments

source
GXBeam.element_curvatureFunction
element_curvature(element, F, M)

Calculate the curvature of a beam element given the resultant force and moments

source
GXBeam.element_propertiesFunction
element_properties(x, icol, beam)
element_properties(x, icol, beam, x0, v0, ω0)
element_properties(x, icol, ibeam, beam, x0, v0, ω0, u, θ, udot, θdot)
element_properties(x, icol, ibeam, beam, x0, v0, ω0, udot, θdot_init, CtCabPdot, CtCabHdot, dt)

Extract/calculate the properties of a specific beam element.

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • x: current state vector
  • icol: starting index for the beam's state variables
  • ibeam: beam element index
  • beam: beam element

Additional Arguments for Dynamic Analyses

  • x0: Global frame origin (for the current time step)
  • v0: Global frame linear velocity (for the current time step)
  • ω0: Global frame angular velocity (for the current time step)

Additional Arguments for Initial Step Analyses

  • u: deflection variables for each beam element
  • θ: rotation variables for each beam element
  • udot: time derivative of u for each beam element
  • θdot: time derivative of θ for each beam element

Additional Arguments for Time Marching Analyses

  • udot_init: 2/dt*u + udot for each beam element from the previous time step
  • θdot_init: 2/dt*θ + θdot for each beam element from the previous time step
  • CtCabPdot_init: 2/dt*C'*Cab*P + C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot_init: 2/dt*C'*Cab*H + C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.element_equationsFunction
element_equations(ΔL, Cab, CtCab, u, θ, F, M, γ, κ)
element_equations(ΔL, Cab, CtCab, u, θ, F, M, γ, κ, v, ω, P, H, V, Ω)
element_equations(ΔL, Cab, CtCab, u, θ, F, M, γ, κ, v, ω, P, H, V, Ω,
    udot, θdot, CtCabPdot, CtCabHdot)

Evaluate the nonlinear equations for a beam element.

There are three implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step or Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments:

  • ΔL: Length of the beam element
  • Ct: Rotation tensor of the beam deformation in the "a" frame, transposed
  • Cab: Direction cosine matrix from "a" to "b" frame for the element
  • CtCab: C'*Cab, precomputed for efficiency
  • u: Displacement variables for the element [u1, u2, u3]
  • θ: Rotation variables for the element [θ1, θ2, θ3]
  • F: Force variables for the element [F1, F2, F3]
  • M: Moment variables for the element [M1, M2, M3]
  • γ: Engineering strains in the element [γ11, 2γ12, 2γ13]
  • κ: Curvatures in the element [κ1, κ2, κ3]

Additional Arguments for Dynamic Analyses

  • v: Linear velocity of element in global frame "a" [v1, v2, v3]
  • ω: Angular velocity of element in global frame "a" [ω1, ω2, ω3]
  • P: Linear momenta for the element [P1, P2, P3]
  • H: Angular momenta for the element [H1, H2, H3]
  • V: Velocity of the element
  • Ω: Rotational velocity of the element

Additional Arguments for Initial Step Analysis

  • udot: user-specified time derivative of u
  • θdot: user-specified time derivative of θ
  • CtCabPdot: C'CabPdot state variable
  • CtCabHdot: C'CabHdot state variable

Additional Arguments for Time Marching Analysis

  • udot: 2/dt*(u-u_p) - udot_p for this beam element where _p denotes values taken from the previous time step
  • θdot: 2/dt*(θ-θ_p) - θdot_p for this beam element where _p denotes values taken from the previous time step
  • CtCabPdot: 2/dt*(C'*Cab*P - (C'*Cab*P)_p) - (C'*Cab*Pdot)_p for this beam element where _p denotes values taken from the previous time step
  • CtCabHdot: 2/dt*(C'*Cab*H - (C'*Cab*H)_p) - (C'*Cab*Hdot)_p for this beam element where _p denotes values taken from the previous time step
source
GXBeam.insert_element_residual!Function
insert_element_residual!(resid, irow_b, irow_b1, irow_p1, irow_b2, irow_p2, f_u1, f_u2, f_ψ1,
    f_ψ2, f_F1, f_F2, f_M1, f_M2)
insert_element_residual!(resid, irow_b, irow_b1, irow_p1, irow_b2, irow_p2, f_u1, f_u2, f_ψ1,
    f_ψ2, f_F1, f_F2, f_M1, f_M2, f_P, f_H)

Insert beam element resultants into the residual equation. Initialize equilibrium and constitutive equations if they are not yet initialized.

If irow_b1 != irow_p1 and/or irow_b2 != irow_p2, assume the equilibrium equations for the left and/or right side are already initialized

There are two implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State, Initial Step, or Time Marching

Arguments

  • resid: System residual vector
  • irow_b1: Row index of the first equation for the left side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p1: Row index of the first equation for the point on the left side of the beam element
  • irow_b1: Row index of the first equation for the right side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p2: Row index of the first equation for the point on the right side of the beam element
  • f_u1, f_u2: Resultant displacements for the left and right side of the beam element, respectively
  • f_ψ1, f_ψ2: Resultant rotations for the left and right side of the beam element, respectively
  • f_F1, f_F2: Resultant forces for the left and right side of the beam element, respectively
  • f_M1, f_M2: Resultant moments for the left and right side of the beam element, respectively

Additional Arguments for Dynamic Analyses

  • f_P: Resultant linear momenta of the beam element
  • f_H: Resultant angular momenta of the beam element
source
GXBeam.element_residual!Function
element_residual!(resid, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2)
element_residual!(resid, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0)
element_residual!(resid, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0, u0, θ0, udot0, θdot0)
element_residual!(resid, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0, udot_init, θdot_init,
    CtCabPdot_init, CtCabHdot_init, dt)

Compute and add a beam element's contributions to the residual vector

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • resid: System residual vector
  • x: current state vector
  • ibeam: beam element index
  • beam: beam element
  • distributed_loads: dictionary with all distributed loads
  • istep: current time step
  • icol: starting index for the beam's state variables
  • irow_b1: Row index of the first equation for the left side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p1: Row index of the first equation for the point on the left side of the beam element
  • irow_b2: Row index of the first equation for the right side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p2: Row index of the first equation for the point on the right side of the beam element

Additional Arguments for Dynamic Analyses

  • x0: Global frame origin (for the current time step)
  • v0: Global frame linear velocity (for the current time step)
  • ω0: Global frame angular velocity (for the current time step)

Additional Arguments for Initial Step Analyses

  • u0: initial deflection variables for each beam element
  • θ0: initial rotation variables for each beam element
  • udot0: initial time derivative of u for each beam element
  • θdot0: initial time derivative of θ for each beam element

Additional Arguments for Time Marching Analyses

  • udot_init: 2/dt*u + udot for each beam element from the previous time step
  • θdot_init: 2/dt*θ + θdot for each beam element from the previous time step
  • CtCabPdot_init: 2/dt*C'*Cab*P + C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot_init: 2/dt*C'*Cab*H + C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.element_jacobian_equationsFunction
element_jacobian_equations(beam, ΔL, Ct, Cab, CtCab, θ, F, M, γ, κ, Ct_θ1,
    Ct_θ2, Ct_θ3)
element_jacobian_equations(beam, ΔL, Cab, CtCab, θ, F, M, γ, κ, ω, P, H, V,
    Ct_θ1, Ct_θ2, Ct_θ3)
element_jacobian_equations(beam, ΔL, Cab, CtCab, θ, F, γ, ω, P, V, θdot)
element_jacobian_equations(beam, ΔL, Cab, CtCab, θ, F, M, γ, κ, ω, P, H, V,
    θdot, dt, Ct_θ1, Ct_θ2, Ct_θ3)

Find the jacobians of the nonlinear equations for a beam element with respect to the state variables given the distributed loads on the beam element and the beam element's properties.

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments:

  • beam: Beam element
  • ΔL: Length of the beam element
  • Ct: Rotation tensor of the beam deformation in the "a" frame, transposed
  • Cab: Direction cosine matrix from "a" to "b" frame for the element
  • CtCab: C'*Cab, precomputed for efficiency
  • u: Displacement variables for the element [u1, u2, u3]
  • θ: Rotation variables for the element [θ1, θ2, θ3]
  • F: Force variables for the element [F1, F2, F3]
  • M: Moment variables for the element [M1, M2, M3]
  • γ: Engineering strains in the element [γ11, 2γ12, 2γ13]
  • κ: Curvatures in the element [κ1, κ2, κ3]
  • Ct_θ1: Gradient of Ct w.r.t. θ[1]
  • Ct_θ2: Gradient of Ct w.r.t. θ[2]
  • Ct_θ3: Gradient of Ct w.r.t. θ[3]
  • f1_θ: Gradient w.r.t. θ of integrated distributed force for left side of element
  • m1_θ: Gradient w.r.t. θ of integrated distributed moment for left side of element
  • f2_θ: Gradient w.r.t. θ of integrated distributed force for right side of element
  • m2_θ: Gradient w.r.t. θ of integrated distributed moment for right side of element

Additional Arguments for Dynamic Analyses

  • v: Linear velocity of element in global frame "a" [v1, v2, v3]
  • ω: Angular velocity of element in global frame "a" [ω1, ω2, ω3]
  • P: Linear momenta for the element [P1, P2, P3]
  • H: Angular momenta for the element [H1, H2, H3]
  • V: Velocity of the element
  • Ω: Rotational velocity of the element

Additional Arguments for Initial Step Analyses

  • udot: user-specified time derivative of u for this beam element
  • θdot: user-specified time derivative of θ for this beam element
  • CtCabPdot: C'*Cab*Pdot (which is a state variable for the initial step analysis)
  • CtCabHdot: C'*Cab*Hdot (which is a state variable for the initial step analysis)

Additional Arguments for Time Marching Analyses

  • udot: -2/dt*u - udot for each beam element from the previous time step
  • θdot_init: -2/dt*θ - θdot for each beam element from the previous time step
  • CtCabPdot: -2/dt*C'*Cab*P - C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot: -2/dt*C'*Cab*H - C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.insert_element_jacobian!Function
insert_element_jacobian!(jacob, irow_b1, irow_p1, irow_b2, irow_p2, icol,
    f_u1_θ, f_u2_θ, f_u1_F, f_u2_F,
    f_ψ1_θ, f_ψ2_θ, f_ψ1_F, f_ψ2_F, f_ψ1_M, f_ψ2_M,
    f_F1_u, f_F2_u, f_F1_θ, f_F2_θ, f_F1_F, f_F2_F, f_F1_M, f_F2_M,
    f_M1_θ, f_M2_θ, f_M1_F, f_M2_F, f_M1_M, f_M2_M)
insert_element_jacobian!(jacob, irow_b1, irow_p1, irow_b2, irow_p2, icol,
    f_u1_θ, f_u2_θ, f_u1_F, f_u2_F, f_u1_P, f_u2_P,
    f_ψ1_θ, f_ψ2_θ, f_ψ1_F, f_ψ2_F, f_ψ1_M, f_ψ2_M, f_ψ1_P, f_ψ2_P, f_ψ1_H, f_ψ2_H,
    f_F1_u, f_F2_u, f_F1_θ, f_F2_θ, f_F1_F, f_F2_F, f_F1_M, f_F2_M,
    f_M1_θ, f_M2_θ, f_M1_F, f_M2_F, f_M1_M, f_M2_M,
    f_P_u, f_P_θ, f_P_P, f_P_H,
    f_H_θ, f_H_P, f_H_H)
insert_element_jacobian!(jacob, irow_b1, irow_p1, irow_b2, irow_p2, icol,
    f_u1_CtCabPdot, f_u2_CtCabPdot, f_u1_F, f_u2_F, f_u1_P, f_u2_P,
    f_ψ1_CtCabHdot, f_ψ2_CtCabHdot, f_ψ1_F, f_ψ2_F, f_ψ1_M, f_ψ2_M, f_ψ1_P, f_ψ2_P, f_ψ1_H, f_ψ2_H,
    f_F1_F, f_F2_F, f_F1_M, f_F2_M,
    f_M1_F, f_M2_F, f_M1_M, f_M2_M,
    f_P_P, f_P_H,
    f_H_P, f_H_H)

Insert the the beam element jacobian entries into the jacobian matrix

There are three implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State or Time Marching
  • Dynamic - Initial Step (for initializing time domain simulations)

Arguments

  • jacob: System jacobian matrix
  • icol_p1: Row/column index of the first unknown for the left endpoint (a value <= 0 indicates the unknowns have been eliminated from the system of equations)
  • irow_b1: Row index of the first equation for the left side of the beam
  • irow_p1: Row index of the first equation for the point on the left side of the beam
  • icol_p2: Row/column index of the first unknown for the right endpoint (a value <= 0 indicates the unknowns have been eliminated from the system of equations)
  • irow_b2: Row index of the first equation for the right side of the beam
  • irow_p2: Row index of the first equation for the point on the right side of the beam
  • icol: Row/Column index corresponding to the first beam state variable

All other arguments use the following naming convention:

  • f_y_x: Jacobian of element equation "y" with respect to state variable "x"
source
GXBeam.element_jacobian!Function
element_jacobian!(jacob, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2)
element_jacobian!(jacob, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0)
element_jacobian!(jacob, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0, u0, θ0, udot0, θdot0)
element_jacobian!(jacob, x, ibeam, beam, distributed_loads, icol, irow_b,
    irow_b1, irow_p1, irow_b2, irow_p2, x0, v0, ω0, udot_init, θdot_init,
    CtCabPdot_init, CtCabHdot_init, dt)

Adds a beam element's contributions to the jacobian matrix

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • jacob: System jacobian matrix
  • x: current state vector
  • ibeam: beam element index
  • beam: beam element
  • distributed_loads: dictionary with all distributed loads
  • istep: current time step
  • icol: starting index for the beam's state variables
  • irow_b1: Row index of the first equation for the left side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p1: Row index of the first equation for the point on the left side of the beam element
  • irow_b2: Row index of the first equation for the right side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p2: Row index of the first equation for the point on the right side of the beam element

Additional Arguments for Dynamic Analyses

  • x0: Global frame origin (for the current time step)
  • v0: Global frame linear velocity (for the current time step)
  • ω0: Global frame angular velocity (for the current time step)

Additional Arguments for Initial Step Analyses

  • u0: initial deflection variables for each beam element
  • θ0: initial rotation variables for each beam element
  • udot0: initial time derivative of u for each beam element
  • θdot0: initial time derivative of θ for each beam element

Additional Arguments for Time Marching Analyses

  • udot_init: 2/dt*u + udot for each beam element from the previous time step
  • θdot_init: 2/dt*θ + θdot for each beam element from the previous time step
  • CtCabPdot_init: 2/dt*C'*Cab*P + C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot_init: 2/dt*C'*Cab*H + C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.element_mass_matrix_propertiesFunction
element_mass_matrix_properties(x, icol, beam)

Extract/Compute the properties needed for mass matrix construction: ΔL, Ct, Cab, CtCab, θ, P, H, Ctdot_cdot1, Ctdot_cdot2, and Ctdot_cdot3

source
GXBeam.element_mass_matrix_equationsFunction
element_mass_matrix_equations(ΔL, Ct, Cab, CtCab, θ, P, H)

Calculates the jacobians of the nonlinear equations for a beam element with respect to the time derivatives of the state variables.

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments:

  • ΔL: Length of the beam element
  • Ct: Rotation tensor of the beam deformation in the "a" frame, transposed
  • Cab: Direction cosine matrix from "a" to "b" frame for the element
  • θ: Rotation variables for the element [θ1, θ2, θ3]
  • P: Linear momenta for the element [P1, P2, P3]
  • H: Angular momenta for the element [H1, H2, H3]
source
GXBeam.insert_element_mass_matrix!Function
insert_element_mass_matrix!(jacob, irow_b1, irow_p1, irow_b2, irow_p2, icol,
    f_u1_θdot, f_u2_θdot, f_u1_Pdot, f_u2_Pdot, f_ψ1_θdot, f_ψ2_θdot,
    f_ψ1_Hdot, f_ψ2_Hdot, f_P_udot, f_H_θdot)

Insert the beam element's contributions into the "mass matrix": the jacobian of the residual equations with respect to the time derivatives of the state variables

Arguments

  • jacob: System jacobian matrix (mass matrix)
  • irow_b1: Row index of the first equation for the left side of the beam
  • irow_p1: Row index of the first equation for the point on the left side of the beam
  • irow_b2: Row index of the first equation for the right side of the beam
  • irow_p2: Row index of the first equation for the point on the right side of the beam
  • icol: Row/Column index corresponding to the first beam state variable

All other arguments use the following naming convention:

  • f_y_x: Jacobian of element equation "y" with respect to state variable "x"
source
GXBeam.element_mass_matrix!Function
element_mass_matrix!(jacob, x, beam, icol, irow_b, irow_b1, irow_p1, irow_b2, irow_p2)

Add the beam element's contributions to the "mass matrix": the jacobian of the residual equations with respect to the time derivatives of the state variables

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • jacob: System jacobian matrix
  • x: current state vector
  • beam: beam element
  • icol: starting index for the beam's state variables
  • irow_b1: Row index of the first equation for the left side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p1: Row index of the first equation for the point on the left side of the beam element
  • irow_b2: Row index of the first equation for the right side of the beam element (a value <= 0 indicates the equations have been eliminated from the system of equations)
  • irow_p2: Row index of the first equation for the point on the right side of the beam element
source
GXBeam.ElementStateType
ElementState

Holds the state variables for an element

Fields:

  • u: Displacement variables for the element (in the global coordinate frame)
  • theta: Wiener-Milenkovic rotational displacement variables for the element
  • F: Resultant forces for the element (in the deformed beam coordinate frame)
  • M: Resultant moments for the element (in the deformed beam coordinate frame)
  • P: Linear momenta of the element (in the deformed beam coordinate frame)
  • H: Angular momenta of the element (in the deformed beam coordinate frame)
source

Loads

GXBeam.PrescribedConditionsType
PrescribedConditions{T}

Describes the forces, moments, displacements, and/or rotations prescribed at a point for each time step

source
GXBeam.DistributedLoadsType
DistributedLoads{T}

Contains the integrated distributed forces and moments for each beam element for each time step.

Fields

  • f1: Integrated non-follower distributed force for the beam element's left endpoint for each time step
  • f2: Integrated non-follower distributed force for the beam element's right endpoint for each time step
  • m1: Integrated non-follower distributed moment for the beam element's left endpoint for each time step
  • m2: Integrated non-follower distributed moment for the beam element's right endpoint for each time step
  • f1_follower: Integrated follower distributed force for the beam element's left endpoint for each time step
  • f2_follower: Integrated follower distributed force for the beam element's right endpoint for each time step
  • m1_follower: Integrated follower distributed moment for the beam element's left endpoint for each time step
  • m2_follower: Integrated follower distributed moment for the beam element's right endpoint for each time step
source

System

GXBeam.AssemblyType
Assembly{TF, TP<:AbstractVector{<:AbstractVector{TF}}, TC<:AbstractVector{<:Integer}, TE<:AbstractVector{Element{TF}}}

Composite type that defines an assembly of connected nonlinear beam elements for analysis.

Fields

  • points: Array of all beam element endpoints
  • start: Array containing point index where each beam element starts
  • stop: Array containing point index where each beam element stops
  • element: Array of Elements
source
GXBeam.curve_triadFunction
curve_triad(Cab, k, s)
curve_triad(Cab, kkt, ktilde, kn, s)

Return the rotation matrix at s along the length of the beam given the curvature vector k and the initial rotation 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 rotation matrix Cab, and curvature vector k.

source
GXBeam.SystemType
System{TF, TV<:AbstractVector{TF}, TM<:AbstractMatrix{TF}, TTF<:AbstractVector{TF}}

Contains the system state, residual vector, and jacobian matrices as well as pointers to be able to access their contents. Also contains additional storage needed for time domain simulations.

Fields:

  • static: Flag indicating whether system matrices are only valid for static analyses
  • x: State vector
  • r: Residual vector
  • K: System jacobian matrix with respect to the state variables
  • M: System jacobian matrix with respect to the time derivative of the state variables
  • irow_pt: Row index of first equilibrium equation for each point
  • irow_beam: Row index of first equation for just this beam element
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
  • icol_pt: Row/Column index of first state variable for each point
  • icol_beam: Row/Column index of first state variable for each beam element
  • current_step: Current time step
  • udot_init: 2/dt*u + udot for each beam element from the previous time step
  • θdot_init: 2/dt*θ + θdot for each beam element from the previous time step
  • CtCabPdot_init: 2/dt*C'*Cab*P + C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot_init: 2/dt*C'*Cab*H + C'*Cab*Hdot for each beam element from the previous time step
source
GXBeam.system_indicesFunction
system_indices(assembly, points, n_connections, static)

Solve for the row indices of the first equilibrium or compatability equations for each point and side of each beam element. Also solve for the row/column index of each point and beam state variable.

Note that this function includes the following logic which reduces the size of the system of equations where possible (without sacrificing any accuracy):

If only two beams meet at a point, the 6 unknowns associated with that point as well as the 6 compatability equations are eliminated from the system, except if specified in the array points. Points for which unknowns have been eliminated are assigned a column index of -1. Beams for which the compatability equations have been eliminated are also assigned an index of -1

Arguments:

  • assembly: Assembly of rigidly connected nonlinear beam elements
  • points: Point indices which should be preserved in the system matrices
  • n_connections: Number of connections to each point
  • static: flag indicating whether analysis is static

Return Arguments:

  • n: total number of equations/unknowns in the system
  • irow_pt: Row index of first equilibrium equation for each point
  • irow_beam: Row index of first equation for just this beam element
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
  • icol_pt: Column index of first state variable for each point
  • icol_beam: Column index of first state variable for each beam element
source
GXBeam.system_residual!Function
system_residual!(resid, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam)
system_residual!(resid, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam, x0, v0, ω0)
system_residual!(resid, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam, x0, v0, ω0,
    u, θ, udot, θdot)
system_residual!(resid, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam, x0, v0, ω0,
    udot, θdot_init, CtCabPdot, CtCabHdot, dt)

Populate the residual vector resid with the results of the residual equations for the system.

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • resid: System residual vector
  • x: Current state variables of the system
  • assembly: Assembly of rigidly connected nonlinear beam elements
  • prescribed_conditions: Dictionary of prescribed conditions at all time steps
  • distributed_loads: Dictionary of distributed loads at all time steps
  • istep: Current time step
  • irow_pt: Row index of first equilibrium equation for each point
  • irow_beam: Row index of first equation for just this beam element
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
  • icol_pt: Column index of first state variable for each point
  • icol_beam: Column index of first state variable for each beam element

Additional Arguments for Dynamic Analyses

  • x0: Global frame origin (for the current time step)
  • v0: Global frame linear velocity (for the current time step)
  • ω0: Global frame angular velocity (for the current time step)

Additional Arguments for Initial Step Analyses

  • u: deflection variables for each beam element
  • θ: rotation variables for each beam element
  • udot: time derivative of u for each beam element
  • θdot: time derivative of θ for each beam element

Additional Arguments for Time Marching Analyses

  • udot: -2/dt*u - udot for each beam element from the previous time step
  • θdot_init: -2/dt*θ - θdot for each beam element from the previous time step
  • CtCabPdot: -2/dt*C'*Cab*P - C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot: -2/dt*C'*Cab*H - C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.system_jacobian!Function
system_jacobian!(jacob, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam)
system_jacobian!(jacob, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam,
    x0, v0, ω0)
system_jacobian!(jacob, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam,
    x0, v0, ω0, u, θ, udot, θdot)
system_jacobian!(jacob, x, assembly, prescribed_conditions, distributed_loads,
    istep, irow_pt, irow_beam, irow_beam1, irow_beam2, icol_pt, icol_beam,
    x0, v0, ω0, udot_init, θdot_init, CtCabPdot_init, CtCabHdot_init, dt)

Populate the jacobian matrix jacob with the jacobian of the residual vector with respect to the state variables.

There are four implementations corresponding to the following analysis types:

  • Static
  • Dynamic - Steady State
  • Dynamic - Initial Step (for initializing time domain simulations)
  • Dynamic - Time Marching

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • jacob: Jacobian matrix
  • x: Vector containing current state variables of the system
  • assembly: Assembly of rigidly connected nonlinear beam elements
  • prescribed_conditions: Dictionary of prescribed conditions at all time steps
  • distributed_loads: Dictionary of distributed loads at all time steps
  • istep: Current time step
  • irow_pt: Row index of first equilibrium equation for each point
  • irow_beam: Row index of first equation for just this beam element
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
  • icol_pt: Column index of first state variable for each point
  • icol_beam: Column index of first state variable for each beam element

Additional Arguments for Dynamic Analyses

  • x0: Global frame origin (for the current time step)
  • v0: Global frame linear velocity (for the current time step)
  • ω0: Global frame angular velocity (for the current time step)

Additional Arguments for Initial Step Analyses

  • u: deflection variables for each beam element
  • θ: rotation variables for each beam element
  • udot: time derivative of u for each beam element
  • θdot: time derivative of θ for each beam element

Additional Arguments for Time Marching Analyses

  • udot: -2/dt*u - udot for each beam element from the previous time step
  • θdot_init: -2/dt*θ - θdot for each beam element from the previous time step
  • CtCabPdot: -2/dt*C'*Cab*P - C'*Cab*Pdot for each beam element from the previous time step
  • CtCabHdot: -2/dt*C'*Cab*H - C'*Cab*Hdot for each beam element from the previous time step
  • dt: time step size
source
GXBeam.system_mass_matrix!Function
system_mass_matrix!(jacob, x, assembly, irow_pt, irow_beam, irow_beam1,
    irow_beam2, icol_pt, icol_beam)

Populate the system "mass matrix", the jacobian of the residual vector with respect to the time derivatives of the state variables.

See "GEBT: A general-purpose nonlinear analysis tool for composite beams" by Wenbin Yu and "Geometrically nonlinear analysis of composite beams using Wiener-Milenković parameters" by Qi Wang and Wenbin Yu.

Arguments

  • jacob: Jacobian matrix
  • x: Vector containing current state variables of the system
  • assembly: Assembly of rigidly connected nonlinear beam elements
  • irow_pt: Row index of first equilibrium equation for each point
  • irow_beam: Row index of first equation for just this beam element
  • irow_beam1: Row index of first equation for the left side of each beam
  • irow_beam2: Row index of first equation for the right side of each beam
  • icol_pt: Column index of first state variable for each point
  • icol_beam: Column index of first state variable for each beam element
source

Index