Library
Public API
Creating an Assembly
GXBeam.curve_length
— Functioncurve_length(start, stop, curvature)
Calculate the length of a curve given its endpoints and its curvature vector
GXBeam.discretize_beam
— Functiondiscretize_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 lengthstart
: Beam starting pointdiscretization
: 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
GXBeam.Assembly
— MethodAssembly(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 endpointsstart
: Array containing point indices where each beam element startsstop
: 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 providingcompliance
compliance
: Array of (6 x 6) compliance matrices for each beam element, defaults tozeros(6,6)
for each beam elementmass
: Array of (6 x 6) mass matrices for each beam element, alternative to providingminv
minv
: Array of (6 x 6) mass matrices for each beam element, defaults to the identity matrix for each beam elementframes
: 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 elementlengths
: Array containing the length of each beam, defaults to the distance between beam endpointsmidpoints
: Array containing the midpoint of each beam element, defaults to the average of the beam element endpoints
Defining Distributed Loads
GXBeam.DistributedLoads
— MethodDistributedLoads(assembly, ibeam[, dt]; kwargs...)
Integrates the specified distributed loads over the element for each time step.
Arguments
assembly
: The beam element assemblyibeam
: The index of the beam element which the distributed load is assigned todt
: Time step size. If omitted a single time step is assumed and specified functions become a function ofs
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 vectormethod = (f, a, b) -> gauss_quadrature(f, a, b)
: Method which integrates functionf
froma
tob
. 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-directionfy = (s, t) -> 0.0
: Distributed non-follower force on beam element in y-directionfz = (s, t) -> 0.0
: Distributed non-follower force on beam element in z-directionmx = (s, t) -> 0.0
: Distributed non-follower moment on beam element in x-directionmy = (s, t) -> 0.0
: Distributed non-follower moment on beam element in y-directionmz = (s, t) -> 0.0
: Distributed non-follower moment on beam element in z-directionfx_follower = (s, t) -> 0.0
: Distributed follower force on beam element in x-directionfy_follower = (s, t) -> 0.0
: Distributed follower force on beam element in y-directionfz_follower = (s, t) -> 0.0
: Distributed follower force on beam element in z-directionmx_follower = (s, t) -> 0.0
: Distributed follower moment on beam element in x-directionmy_follower = (s, t) -> 0.0
: Distributed follower moment on beam element in y-directionmz_follower = (s, t) -> 0.0
: Distributed follower moment on beam element in z-direction
Defining Prescribed Conditions
GXBeam.PrescribedConditions
— MethodPrescribedConditions(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 vectorux
: Prescribed x-direction displacement of the pointuy
: Prescribed y-direction displacement of the pointuz
: Prescribed z-direction displacement of the pointtheta_x
: Prescribed first Wiener-Milenkovic parameter of the pointtheta_y
: Prescribed second Wiener-Milenkovic parameter of the pointtheta_z
: Prescribed third Wiener-Milenkovic parameter of the pointFx
: Prescribed force in x-direction applied on the pointFy
: Prescribed force in y-direction applied on the pointFz
: Prescribed force in z-direction applied on the pointMx
: Prescribed moment about x-axis applied on the pointMy
: Prescribed moment about y-axis applied on the pointMz
: Prescribed moment about z-axis applied on the pointFx_follower
: Prescribed follower force in x-direction applied on the pointFy_follower
: Prescribed follower force in y-direction applied on the pointFz_follower
: Prescribed follower force in z-direction applied on the pointMx_follower
: Prescribed follower moment about x-axis applied on the pointMy_follower
: Prescribed follower moment about y-axis applied on the pointMz_follower
: Prescribed follower moment about z-axis applied on the point
Pre-Initializing Memory for an Analysis
GXBeam.System
— MethodSystem([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 resultingSystem
objectassembly
: Assembly of rigidly connected nonlinear beam elementspoints
: 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
GXBeam.reset_state!
— Functionreset_state!(system)
Reset the state variables in system
(stored in system.x
) to zero.
Performing an Analysis
GXBeam.static_analysis
— Functionstatic_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 inkeys(prescribed_conditions)
distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: Dictionary holding DistributedLoads composite types for the beam elements inkeys(distributed_loads)
linear = false
: Set totrue
for a linear analysismethod = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve nonlinear system of equationsftol = 1e-9
: tolerance for solving nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equationsnstep = 1
: Number of time steps. May be used in conjunction with time varying prescribed conditions and distributed loads to gradually increase displacements/loads.
GXBeam.static_analysis!
— Functionstatic_analysis!(system, assembly; kwargs...)
Pre-allocated version of static_analysis
. Uses the state variables stored in system
as an initial guess for iterating.
GXBeam.steady_state_analysis
— Functionsteady_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 inkeys(prescribed_conditions)
distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: Dictionary holding DistributedLoads composite types for the beam elements inkeys(distributed_loads)
linear = false
: Set totrue
for a linear analysismethod = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve nonlinear system of equationsftol = 1e-9
: tolerance for solving nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equationsnstep = 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 originlinear_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
GXBeam.steady_state_analysis!
— Functionsteady_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.
GXBeam.eigenvalue_analysis
— Functioneigenvalue_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 inkeys(prescribed_conditions)
distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: Dictionary holding DistributedLoads composite types for the beam elements inkeys(distributed_loads)
linear = false
: Set totrue
for a linear analysismethod = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve nonlinear system of equationsftol = 1e-9
: tolerance for solving nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equationsnstep = 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 originlinear_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 stepnev = 6
: Number of eigenvalues to compute
GXBeam.eigenvalue_analysis!
— Functioneigenvalue_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.
GXBeam.time_domain_analysis
— Functiontime_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 inkeys(prescribed_conditions)
distributed_loads = Dict{Int,DistributedLoads{Float64}}()
: Dictionary holding DistributedLoads composite types for the beam elements inkeys(distributed_loads)
linear = false
: Set totrue
for a linear analysismethod = :newton
: Method (as defined in NLsolve) to solve nonlinear system of equationslinesearch = LineSearches.LineSearches.BackTracking(maxstep=1e6)
: Line search used to solve nonlinear system of equationsftol = 1e-9
: tolerance for solving nonlinear system of equationsiterations = 1000
: maximum iterations for solving the nonlinear system of equationsnstep = 1
: The total length of the time vectororigin = zeros(3)
: Global frame originlinear_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 tou
thetadot0=fill(zeros(3), length(assembly.elements))
: initial time derivative with respect totheta
save=1:nstep
: Steps at which to save the time history
GXBeam.time_domain_analysis!
— Functiontime_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.
Post-Processing
GXBeam.AssemblyState
— MethodAssemblyState(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.
GXBeam.AssemblyState
— TypeAssemblyState{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 ofPointState
s for each point in the assemblyelements::TE
: Array ofElementState
s for each element in the assembly
GXBeam.left_eigenvectors
— Functionleft_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.
GXBeam.correlate_eigenmodes
— Functioncorrelate_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.
GXBeam.wiener_milenkovic
— Functionwiener_milenkovic(c)
Construct a Wiener-Milenkovic rotation matrix, given the three Wiener-Milenkovic parameters in c
.
GXBeam.write_vtk
— Functionwrite_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)
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"
Private API
Math
GXBeam.tilde
— Functiontilde(x)
Construct the cross product operator matrix
GXBeam.rotation_parameter_scaling
— Functionrotation_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.
GXBeam.get_C
— Functionget_C(θ)
Returns the rotation matrix C
given the three angular parameters in θ
.
GXBeam.get_C_t
— Functionget_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
.
GXBeam.get_C_θ
— Functionget_C_θ([C, ] θ)
Calculate the derivative of the Wiener-Milenkovic rotation matrix C
with respect to each of the rotation parameters in θ
.
GXBeam.get_C_θdot
— Functionget_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.
GXBeam.get_Q
— Functionget_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 θ
.
GXBeam.get_Q_θ
— Functionget_Q_θ(θ)
get_Q_θ(Q, θ)
Calculate the derivative of the matrix Q
with respect to each of the rotation parameters in c
.
GXBeam.get_Qinv
— Functionget_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 θ
.
GXBeam.get_Qinv_θ
— Functionget_Qinv_θ(θ)
Calculate the derivative of the matrix inverse Qinv
with respect to each of the rotation parameters in θ
.
GXBeam.mul3
— Functionmul3(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
.
GXBeam.gauss_quadrature
— Functiongauss_quadrature(f, a, b)
Default gauss-quadrature function used for integrating distributed loads.
Points
GXBeam.point_variables
— Functionpoint_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])
GXBeam.insert_point_residual!
— Functioninsert_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))
GXBeam.point_residual!
— Functionpoint_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 vectorx
: current state vectoripoint
: index of pointassembly
: assembly of interconnected beam elementsprescribed_conditions
: dictionary of prescribed conditionsicol
: starting index for the point's state variablesirow_p
: Row index of the first equilibrium equation for the pointirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beam
GXBeam.point_follower_jacobians
— Functionpoint_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
GXBeam.insert_point_jacobian!
— Functioninsert_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))
GXBeam.point_jacobian!
— Functionpoint_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 vectorx
: current state vectoripoint
: index of pointassembly
: assembly of interconnected beam elementsprescribed_conditions
: dictionary of prescribed conditionsistep
: current time stepicol
: starting index for the point's state variablesirow_p
: Row index of the first equilibrium equation for the pointirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beam
GXBeam.PointState
— TypePointState
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 pointF
: Externally applied forces on the point (in the global coordinate frame)M
: Externally applied moments on the point (in the global coordinate frame)
Elements
GXBeam.Element
— TypeElement{TF}
Composite type that defines a beam element's properties
Fields
L
: Length of the beam elementx
: Location of the beam element (the midpoint of the beam element)C11
: Upper left portion of the beam element's compliance matrixC12
: Upper right portion of the beam element's compliance matrixC22
: Lower right portion of the beam element's compliance matrixminv11
: Upper left portion of the inverse of the beam element's mass matrixminv12
: Upper right portion of the inverse of the beam element's mass matrixminv22
: Lower right portion of the inverse of the beam element's mass matrixCab
: Rotation matrix from the global frame to beam element's frame
GXBeam.element_strain
— Functionelement_strain(element, F, M)
Calculate the strain of a beam element given the resultant force and moments
GXBeam.element_curvature
— Functionelement_curvature(element, F, M)
Calculate the curvature of a beam element given the resultant force and moments
GXBeam.element_linear_velocity
— Functionelement_linear_velocity(element, P, H)
Calculate the linear velocity of a beam element given the linear and angular momenta
GXBeam.element_angular_velocity
— Functionelement_angular_velocity(element, P, H)
Calculate the angular velocity of a beam element given the linear and angular momenta
GXBeam.element_properties
— Functionelement_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 vectoricol
: starting index for the beam's state variablesibeam
: beam element indexbeam
: 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 elementudot
: 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 stepCtCabPdot_init
:2/dt*C'*Cab*P + C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot_init
:2/dt*C'*Cab*H + C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.dynamic_element_properties
— Functiondynamic_element_properties(x, icol, beam, x0, v0, ω0)
Extract/Compute v
, ω
, P
, H
, V
, and Ω
.
GXBeam.element_equations
— Functionelement_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 elementCt
: Rotation tensor of the beam deformation in the "a" frame, transposedCab
: Direction cosine matrix from "a" to "b" frame for the elementCtCab
:C'*Cab
, precomputed for efficiencyu
: 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 variableCtCabHdot
: 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 stepCtCabPdot
: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 stepCtCabHdot
: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
GXBeam.insert_element_residual!
— Functioninsert_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 vectorirow_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 elementirow_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 elementf_u1
,f_u2
: Resultant displacements for the left and right side of the beam element, respectivelyf_ψ1
,f_ψ2
: Resultant rotations for the left and right side of the beam element, respectivelyf_F1
,f_F2
: Resultant forces for the left and right side of the beam element, respectivelyf_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 elementf_H
: Resultant angular momenta of the beam element
GXBeam.element_residual!
— Functionelement_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 vectorx
: current state vectoribeam
: beam element indexbeam
: beam elementdistributed_loads
: dictionary with all distributed loadsistep
: current time stepicol
: starting index for the beam's state variablesirow_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 elementirow_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 elementudot0
: 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 stepCtCabPdot_init
:2/dt*C'*Cab*P + C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot_init
:2/dt*C'*Cab*H + C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.element_jacobian_equations
— Functionelement_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 elementCt
: Rotation tensor of the beam deformation in the "a" frame, transposedCab
: Direction cosine matrix from "a" to "b" frame for the elementCtCab
:C'*Cab
, precomputed for efficiencyu
: 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 ofCt
w.r.t.θ[1]
Ct_θ2
: Gradient ofCt
w.r.t.θ[2]
Ct_θ3
: Gradient ofCt
w.r.t.θ[3]
f1_θ
: Gradient w.r.t. θ of integrated distributed force for left side of elementm1_θ
: Gradient w.r.t. θ of integrated distributed moment for left side of elementf2_θ
: Gradient w.r.t. θ of integrated distributed force for right side of elementm2_θ
: 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 elementCtCabPdot
: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 stepCtCabPdot
:-2/dt*C'*Cab*P - C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot
:-2/dt*C'*Cab*H - C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.insert_element_jacobian!
— Functioninsert_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 matrixicol_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 beamirow_p1
: Row index of the first equation for the point on the left side of the beamicol_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 beamirow_p2
: Row index of the first equation for the point on the right side of the beamicol
: 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"
GXBeam.element_jacobian!
— Functionelement_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 matrixx
: current state vectoribeam
: beam element indexbeam
: beam elementdistributed_loads
: dictionary with all distributed loadsistep
: current time stepicol
: starting index for the beam's state variablesirow_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 elementirow_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 elementudot0
: 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 stepCtCabPdot_init
:2/dt*C'*Cab*P + C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot_init
:2/dt*C'*Cab*H + C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.element_mass_matrix_properties
— Functionelement_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
GXBeam.element_mass_matrix_equations
— Functionelement_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 elementCt
: Rotation tensor of the beam deformation in the "a" frame, transposedCab
: 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]
GXBeam.insert_element_mass_matrix!
— Functioninsert_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 beamirow_p1
: Row index of the first equation for the point on the left side of the beamirow_b2
: Row index of the first equation for the right side of the beamirow_p2
: Row index of the first equation for the point on the right side of the beamicol
: 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"
GXBeam.element_mass_matrix!
— Functionelement_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 matrixx
: current state vectorbeam
: beam elementicol
: starting index for the beam's state variablesirow_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 elementirow_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
GXBeam.ElementState
— TypeElementState
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)
Loads
GXBeam.PrescribedConditions
— TypePrescribedConditions{T}
Describes the forces, moments, displacements, and/or rotations prescribed at a point for each time step
GXBeam.DistributedLoads
— TypeDistributedLoads{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
System
GXBeam.Assembly
— TypeAssembly{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 endpointsstart
: Array containing point index where each beam element startsstop
: Array containing point index where each beam element stopselement
: Array ofElement
s
GXBeam.curve_triad
— Functioncurve_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
.
GXBeam.curve_coordinates
— Functioncurve_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
.
GXBeam.System
— TypeSystem{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 analysesx
: State vectorr
: Residual vectorK
: System jacobian matrix with respect to the state variablesM
: System jacobian matrix with respect to the time derivative of the state variablesirow_pt
: Row index of first equilibrium equation for each pointirow_beam
: Row index of first equation for just this beam elementirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beamicol_pt
: Row/Column index of first state variable for each pointicol_beam
: Row/Column index of first state variable for each beam elementcurrent_step
: Current time stepudot_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 stepCtCabPdot_init
:2/dt*C'*Cab*P + C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot_init
:2/dt*C'*Cab*H + C'*Cab*Hdot
for each beam element from the previous time step
GXBeam.point_connections
— Functionpoint_connections(assembly)
Count the number of beams connected to each point
GXBeam.system_indices
— Functionsystem_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 elementspoints
: Point indices which should be preserved in the system matricesn_connections
: Number of connections to each pointstatic
: flag indicating whether analysis is static
Return Arguments:
n
: total number of equations/unknowns in the systemirow_pt
: Row index of first equilibrium equation for each pointirow_beam
: Row index of first equation for just this beam elementirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beamicol_pt
: Column index of first state variable for each pointicol_beam
: Column index of first state variable for each beam element
GXBeam.system_residual!
— Functionsystem_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 vectorx
: Current state variables of the systemassembly
: Assembly of rigidly connected nonlinear beam elementsprescribed_conditions
: Dictionary of prescribed conditions at all time stepsdistributed_loads
: Dictionary of distributed loads at all time stepsistep
: Current time stepirow_pt
: Row index of first equilibrium equation for each pointirow_beam
: Row index of first equation for just this beam elementirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beamicol_pt
: Column index of first state variable for each pointicol_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 elementudot
: 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 stepCtCabPdot
:-2/dt*C'*Cab*P - C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot
:-2/dt*C'*Cab*H - C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.system_jacobian!
— Functionsystem_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 matrixx
: Vector containing current state variables of the systemassembly
: Assembly of rigidly connected nonlinear beam elementsprescribed_conditions
: Dictionary of prescribed conditions at all time stepsdistributed_loads
: Dictionary of distributed loads at all time stepsistep
: Current time stepirow_pt
: Row index of first equilibrium equation for each pointirow_beam
: Row index of first equation for just this beam elementirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beamicol_pt
: Column index of first state variable for each pointicol_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 elementudot
: 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 stepCtCabPdot
:-2/dt*C'*Cab*P - C'*Cab*Pdot
for each beam element from the previous time stepCtCabHdot
:-2/dt*C'*Cab*H - C'*Cab*Hdot
for each beam element from the previous time stepdt
: time step size
GXBeam.system_mass_matrix!
— Functionsystem_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 matrixx
: Vector containing current state variables of the systemassembly
: Assembly of rigidly connected nonlinear beam elementsirow_pt
: Row index of first equilibrium equation for each pointirow_beam
: Row index of first equation for just this beam elementirow_beam1
: Row index of first equation for the left side of each beamirow_beam2
: Row index of first equation for the right side of each beamicol_pt
: Column index of first state variable for each pointicol_beam
: Column index of first state variable for each beam element
Index
GXBeam.Assembly
GXBeam.Assembly
GXBeam.AssemblyState
GXBeam.AssemblyState
GXBeam.DistributedLoads
GXBeam.DistributedLoads
GXBeam.Element
GXBeam.ElementState
GXBeam.PointState
GXBeam.PrescribedConditions
GXBeam.PrescribedConditions
GXBeam.System
GXBeam.System
GXBeam.correlate_eigenmodes
GXBeam.curve_coordinates
GXBeam.curve_length
GXBeam.curve_triad
GXBeam.discretize_beam
GXBeam.dynamic_element_properties
GXBeam.eigenvalue_analysis
GXBeam.eigenvalue_analysis!
GXBeam.element_angular_velocity
GXBeam.element_curvature
GXBeam.element_equations
GXBeam.element_jacobian!
GXBeam.element_jacobian_equations
GXBeam.element_linear_velocity
GXBeam.element_mass_matrix!
GXBeam.element_mass_matrix_equations
GXBeam.element_mass_matrix_properties
GXBeam.element_properties
GXBeam.element_residual!
GXBeam.element_strain
GXBeam.gauss_quadrature
GXBeam.get_C
GXBeam.get_C_t
GXBeam.get_C_θ
GXBeam.get_C_θdot
GXBeam.get_Q
GXBeam.get_Q_θ
GXBeam.get_Qinv
GXBeam.get_Qinv_θ
GXBeam.insert_element_jacobian!
GXBeam.insert_element_mass_matrix!
GXBeam.insert_element_residual!
GXBeam.insert_point_jacobian!
GXBeam.insert_point_residual!
GXBeam.left_eigenvectors
GXBeam.mul3
GXBeam.point_connections
GXBeam.point_follower_jacobians
GXBeam.point_jacobian!
GXBeam.point_residual!
GXBeam.point_variables
GXBeam.reset_state!
GXBeam.rotation_parameter_scaling
GXBeam.static_analysis
GXBeam.static_analysis!
GXBeam.steady_state_analysis
GXBeam.steady_state_analysis!
GXBeam.system_indices
GXBeam.system_jacobian!
GXBeam.system_mass_matrix!
GXBeam.system_residual!
GXBeam.tilde
GXBeam.time_domain_analysis
GXBeam.time_domain_analysis!
GXBeam.wiener_milenkovic
GXBeam.write_vtk