Similarly to the previous section, the aim is to describe the direct time solvers available for solving Mdof systems. As a remainder, the equation of motion of the system is given by: \[
\mathbf{M} \ddot{\mathbf{x}}(t) + \mathbf{C} \dot{\mathbf{x}}(t) + \mathbf{K} \mathbf{x}(t) = \mathbf{F}(t)
\] where \(\mathbf{M}\) is the mass matrix, \(\mathbf{C}\) is the damping matrix, \(\mathbf{K}\) is the stiffness matrix, \(\mathbf{x}(t)\) is the displacement vector, and \(\mathbf{F}(t)\) is the vector of external forces.
Direct-time integration solvers are used to solve the previous equation in the time domain. In Julia , we have the option of using DifferentialEquations.jl which provides a wide range of solvers for solving differential equations.
The StructuralVibration package provides its own set of fixed-step linear second-order ODE solvers for pedagogical purposes and for the sake of completeness of the package. It should also be noted that some of the solvers implemented in the package are not available in DiffentialEquations.jl. This is the case of the solvers of the Newmark family, which are widely used in structural dynamics.
The solvers available in the StructuralVibration package are:
sol: Solution structure containing the response of the system at the given time points
u: Displacement
du: Velocity
ddu: Acceleration
2 Free response
The free response of a 2-dofs system is obtained by setting the external forces to zero. In this condition, the motion of the system is initiated by non-zero initial conditions. In the following example, the response provided by the direct time solvers is compared with the solution obtained by the modal solver described in Modal time solvers - Section 1 for the purpose of cross-validation.
# System parametersM =Diagonal([2., 1.])K = [6. -2.; -2. 4.]ξ =0.05ω, Φ =eigenmode(K, M)C =modal_damping_matrix(M, ω, ξ, Φ)# Time vectort =0.:1e-2:30.# Initial conditionsx0 = [0.2, 0.1]v0 =zeros(2)u0 = (x0, v0)# External forcesF_free =zeros(2, length(t))# Direct time problemprob_free =DirectTimeProblem(K, M, C, F_free, u0, t)x_free_gα =solve(prob_free).ux_free_cd =solve(prob_free, CentralDiff()).ux_free_rk =solve(prob_free, RK4()).u# Modal time problemprob_free_modal =FreeModalTimeProblem(K, M, ξ, u0, t)x_free_modal =solve(prob_free_modal).u
3 Forced response
The forced response of a 2-dofs system is obtained by applying a haversine excitation to the first dof of the system. As in the previous section, the response provided by the direct time solvers is compared with the solution obtained by the modal solver described in Modal time solvers - Section 2 for the purpose of cross-validation.
# System parametersM =Diagonal([2., 1.])K = [6. -2.; -2. 4.]ξ =0.05ω, Φ =eigenmode(K, M)C =modal_damping_matrix(M, ω, ξ, Φ)# Time vectort =0.:1e-2:30.# Initial conditionsx0 =zeros(2)v0 =zeros(2)u0 = (x0, v0)# External forcesF0 =10.tstart =2.duration =5.haversine =HaverSine(F0, tstart, duration)F0 =excitation(haversine, t)F =zeros(2, length(t))F[1, :] .= F0# Direct time problemprob_forced =DirectTimeProblem(K, M, C, F, u0, t)x_forced_gα =solve(prob_forced).ux_forced_cd =solve(prob_forced, CentralDiff()).ux_forced_rk =solve(prob_forced, RK4()).u# Modal time problemprob_forced_modal =ForcedModalTimeProblem(K, M, ξ, F, u0, t)x_forced_modal =solve(prob_forced_modal).u
Footnotes
N. M. Newmark. “A method of computation for structural dynamics”. Journal of the Engineering Mechanics Division, 85(3), 67-94, 1959.↩︎
H. M. Hilber, T. J. Hughes and R. L. Taylor. “Improved numerical dissipation for time integration algorithms in structural dynamics”. Earthquake Engineering & Structural Dynamics 5(3), 283–292. 1977.↩︎
W. L. Wood, M. Bossak, O. Zienkiewicz. “An alpha modification of Newmark’s method”. International Journal for Numerical Methods in Engineering, 15(10), 1562-1566, 1980.↩︎
L. Fox, E. T. Goodwin. “Some new methods for numerical integration of ordinary differential equations”. Mathematical Proceedings of the Cambridge Philosophical Society, 45(3), 373-388, 1949.↩︎
M. Géradin and D. Rixen. “Mechanical Vibrations: Theory and Application to Structural Dynamics”. Third Edition. John Wiley & Sons, 2015.↩︎
J. Chung and G. M. Hulbert. “A time integration algorithm for structural dynamics with improved numerical dissipation: the generalized-alpha method”. Journal of Applied Mechanics, 60(2), 371-375, 1993.↩︎