Hi,
Got another issue... it concerns problems with "free-flight" boundary conditions, say for a beam free to move in 3D space.
First, I verified what seems to be an issue even in a 2D constrained flight. The problem is the "flying spaghetti" from 1, here is the GXBeam input code:
using GXBeam, LinearAlgebra
# SIMO & VU-QUOC's 2D flying spaghetti
# beam properties
L = 10;
EA = 1e4;
GA = EA;
EI = 500;
GJ = EI;
rhoA = 1;
rhoIs = 20;
rhoI = 10;
# Number of elements
nelem = 20
# Frame for initial rotation of the beam
Cab = [0.6 0.0 0.8
0.0 1.0 0.0
-0.8 0.0 0.6]
# Discretization variables
lengths, points, midpoints, Cab = discretize_beam(L, [0, 0, 0], nelem, frame=Cab)
# index of endpoints of each beam element
start = 1:nelem
stop = 2:nelem+1
# stiffness matrix for each beam element
stiffness = fill(Diagonal([EA, GA, GA, GJ, EI, EI]), nelem)
# mass matrix for each beam element
mass = fill(Diagonal([rhoA, rhoA, rhoA, rhoIs, rhoI, rhoI]), nelem)
# create assembly
assembly = Assembly(points, start, stop;
stiffness = stiffness,
mass = mass,
frames = Cab,
lengths = lengths,
midpoints = midpoints)
# external forces
tau = 2.5
M = (t) -> begin
if 0.0 <= t <= tau
80
else
zero(t)
end
end
# prescribed conditions
prescribed_conditions = (t) -> begin
Dict(
# left side: no out of plane movement
1 => PrescribedConditions(uy=0, theta_z=0),
# right side: no out of plane movement, applied forces
nelem+1 => PrescribedConditions(Fx=M(t)/10, My=M(t), uy=0, theta_z=0)
)
end
# simulation time
t = 0:1e-2:13
# analysis
system, history, converged = time_domain_analysis(assembly, t;
prescribed_conditions = prescribed_conditions)
# Get displacement outputs at both ends of the beam
point = [1 nelem+1]
field = [:u, :u]
direction = [1, 3]
ylabel = ["\$u\$ (\$m\$)","\$w\$ (\$m\$)"]
u1 = [getproperty(state.points[point[1]], field[1])[direction[1]] for state in history]
w1 = [getproperty(state.points[point[1]], field[2])[direction[2]] for state in history]
uend = [getproperty(state.points[point[2]], field[1])[direction[1]] for state in history]
wend = [getproperty(state.points[point[2]], field[2])[direction[2]] for state in history]
# plot results
using Plots
pyplot()
# plot horizontal displacement
plot(xlabel = "Time (s)",
ylabel = ylabel[1],
grid = false,
overwrite_figure=false
)
plot!(t, u1, label="Node 1")
plot!(t, uend, label="Node end")
plot!(show=true)
# plot vertical displacement
plot(xlabel = "Time (s)",
ylabel = ylabel[2],
grid = false,
overwrite_figure=false
)
plot!(t, w1, label="Node 1")
plot!(t, wend, label="Node end")
plot!(show=true)
A motion plot of the problem is as follows:

As you can see, after about 7 seconds, the beam starts to exhibit an unexpected motion, moving towards the left and downwards, even though its momentum should make it keep going to the right indefinitely, with no vertical movement. The authors in the original paper1 do not make it clear if they modified their motion plots in any way, but a more recent work2 confirms that indeed the beam keeps moving to the right:

Setting the prescribed condition u_z = 0 for the beam's midpoint fixes the problem, but such boundary condition is not stated by either 1 or 2, neither should it be necessary. Comparing solutions, the deformed shapes of the beam seem to agree, but the positions in space do not.
I also investigated the 3D version of that problem3, in which an out-of-plane moment is also applied, a true free-flight case. Here is the GXBeam file:
using GXBeam, LinearAlgebra
# SIMO & VU-QUOC's 3D flying spaghetti
# beam properties
L = 10;
EA = 1e4;
GA = EA;
EI = 500;
GJ = EI;
rhoA = 1;
rhoIs = 20;
rhoI = 10;
# Number of elements
nelem = 20
# Frame for initial rotation of the beam
Cab = [0.6 0.0 0.8
0.0 1.0 0.0
-0.8 0.0 0.6]
# Discretization variables
lengths, points, midpoints, Cab = discretize_beam(L, [0, 0, 0], nelem, frame=Cab)
# index of endpoints of each beam element
start = 1:nelem
stop = 2:nelem+1
# stiffness matrix for each beam element
stiffness = fill(Diagonal([EA, GA, GA, GJ, EI, EI]), nelem)
# mass matrix for each beam element
mass = fill(Diagonal([rhoA, rhoA, rhoA, rhoIs, rhoI, rhoI]), nelem)
# create assembly
assembly = Assembly(points, start, stop;
stiffness = stiffness,
mass = mass,
frames = Cab,
lengths = lengths,
midpoints = midpoints)
# forces
M0 = 200
tau = 2.5
M = (t) -> begin
if 0.0 <= t <= tau
M0*t/tau
elseif tau < t <= 2*tau
2*M0*(1.0-t/(2*tau))
else
zero(t)
end
end
# prescribed conditions
prescribed_conditions = (t) -> begin
Dict(
# forces on right side
nelem+1 => PrescribedConditions(Fx=M(t)/10, My=M(t), Mz=M(t)/2)
)
end
# simulation time (problems get serious after about 3 seconds)
t = 0:1e-3:3.14
# analysis
system, history, converged = time_domain_analysis(assembly, t;
prescribed_conditions = prescribed_conditions)
# Get displacement outputs at both ends of the beam
point = [1 nelem+1]
field = [:u, :u, :theta]
direction = [1, 3, 2]
ylabel = ["\$u\$ (\$m\$)","\$w\$ (\$m\$)","Rodrigues parameter \$\\theta_y\$"]
u1 = [getproperty(state.points[point[1]], field[1])[direction[1]] for state in history]
w1 = [getproperty(state.points[point[1]], field[2])[direction[2]] for state in history]
uend = [getproperty(state.points[point[2]], field[1])[direction[1]] for state in history]
wend = [getproperty(state.points[point[2]], field[2])[direction[2]] for state in history]
thetay_end = [getproperty(state.points[point[2]], field[3])[direction[3]] for state in history]
# plot results
using Plots
pyplot()
# plot u
plot(xlabel = "Time (s)",
ylabel = ylabel[1],
grid = true,
overwrite_figure=false
)
plot!(t, u1, label="Node 1")
plot!(t, uend, label="Node end")
plot!(show=true)
# plot w
plot(xlabel = "Time (s)",
ylabel = ylabel[2],
grid = true,
overwrite_figure=false
)
plot!(t, w1, label="Node 1")
plot!(t, wend, label="Node end")
plot!(show=true)
# plot theta_y
plot(xlabel = "Time (s)",
ylabel = ylabel[3],
grid = true,
overwrite_figure=false
)
plot!(t, thetay_end, label="")
plot!(show=true)
The plot of the theta_y rotation parameter reveals an unstable behavior of the beam, and the solution only converges up until around 3.16 seconds:

No discretization or time step refinement could remedy the numerical instability. A plane view motion plot makes it clear that there is an issue in the simulation:

Since the program appears to work perfectly in cases where the structure is restrained in some way, I was led to believe the issue only happens in these free-flight cases. Any ideas of what may be going on? Maybe there is a need to integrate the "a" frame states over time as well?
Hi,
Got another issue... it concerns problems with "free-flight" boundary conditions, say for a beam free to move in 3D space.
First, I verified what seems to be an issue even in a 2D constrained flight. The problem is the "flying spaghetti" from 1, here is the GXBeam input code:
A motion plot of the problem is as follows:
As you can see, after about 7 seconds, the beam starts to exhibit an unexpected motion, moving towards the left and downwards, even though its momentum should make it keep going to the right indefinitely, with no vertical movement. The authors in the original paper1 do not make it clear if they modified their motion plots in any way, but a more recent work2 confirms that indeed the beam keeps moving to the right:
Setting the prescribed condition u_z = 0 for the beam's midpoint fixes the problem, but such boundary condition is not stated by either 1 or 2, neither should it be necessary. Comparing solutions, the deformed shapes of the beam seem to agree, but the positions in space do not.
I also investigated the 3D version of that problem3, in which an out-of-plane moment is also applied, a true free-flight case. Here is the GXBeam file:
The plot of the theta_y rotation parameter reveals an unstable behavior of the beam, and the solution only converges up until around 3.16 seconds:
No discretization or time step refinement could remedy the numerical instability. A plane view motion plot makes it clear that there is an issue in the simulation:
Since the program appears to work perfectly in cases where the structure is restrained in some way, I was led to believe the issue only happens in these free-flight cases. Any ideas of what may be going on? Maybe there is a need to integrate the "a" frame states over time as well?
Footnotes
Simo, J. C., & Vu-Quoc, L. (1986). On the dynamics of flexible beams under large overall motions—The plane case: Part II. Journal of Applied Mechanics, 53, 855-863. ↩ ↩2 ↩3
Palacios, R., & Cea, A. (2019). Nonlinear modal condensation of large finite element models: Application of Hodges’s intrinsic theory. AIAA Journal, 57(10), 4255-4268. ↩ ↩2
Simo, J. C., & Vu-Quoc, L. (1988). On the dynamics in space of rods undergoing large motions—a geometrically exact approach. Computer Methods in Applied Mechanics and Engineering, 66(2), 125-161. ↩