For solving time-domain problems, the StructuralVibration.jl package provides two different approach. The first one consist in solving the problem using the discrete-time state equation (see Mechanical solvers - Section 3.2 for details). In this case, four discretization methods are available: Zero-Order Hold (:zoh), First-Order Hold (:foh), Band-Limited Hold (:blh), and Fourth order Runge-Kutta scheme (:rk4).
The second approach is to solve the problem using the continuous-time state equation. In this case, the user can build a continuous state-space model and solve it using DifferentialEquations.jl to take advantage of one of its numerous solvers or use the Fourth order Runge-Kutta scheme implemented in the package.
The frequency-domain solvers are based on the state-space representation of the system, which includes both the state equation and the output equation. It is given that: \[
\begin{cases}
\dot{\mathbf{z}}(t) &= \mathbf{A}_c \mathbf{z}(t) + \mathbf{B}_c \mathbf{u}(t) \\
\mathbf{y}(t) &= \mathbf{C}_c \mathbf{z}(t) + \mathbf{D}_c \mathbf{u}(t)
\end{cases},
\] where:
\(\mathbf{z}(t)\): State vector
\(\mathbf{u}(t)\): Input vector
\(\mathbf{y}(t)\): Output vector
\(\mathbf{A}_c\): System matrix
\(\mathbf{B}_c\): Input matrix
\(\mathbf{C}_c\): Observation matrix
\(\mathbf{D}_c\): Direct feedthrough matrix
It is also possible to express the state-space system. To this end, the state-space representation is transformed into modal coordinates, i.e., \(\mathbf{z} = \mathbf{\Psi} \mathbf{q}\), where \(\mathbf{\Psi}\) is the matrix of eigenvectors of the system matrix \(\mathbf{A}_c\). The modal state-space representation is given by: \[
\begin{cases}
\dot{\mathbf{q}}(t) &= \mathbf{\Lambda} \mathbf{q}(t) + \mathbf{B}_m \mathbf{u}(t) \\
\mathbf{y}(t) &= \mathbf{C}_m \mathbf{q}(t) + \mathbf{D}_m \mathbf{u}(t)
\end{cases},
\] where:
\(\mathbf{D}_m = \mathbf{D}_c\): Modal direct feedthrough matrix
Finally
2.1 Frequency Response Function
The Frequency Response Function (FRF) is a complex function that describes the steady-state response of a system to a sinusoidal input. The FRF is defined as the ratio of the output to the input in the frequency domain. From the state-space representation, the transfer function of the system at a given angular frequency \(\omega\) can be obtained in a straightforward manner and is given by: \[
\mathbf{H}(\omega) = \mathbf{C}_c (j\omega \mathbf{I} - \mathbf{A}_c)^{-1} \mathbf{B}_c + \mathbf{D}_c,
\] or equivalently from the modes of the system: \[
\mathbf{H}(\omega) = \mathbf{C}_m (j\omega \mathbf{I} - \mathbf{\Lambda} )^{-1} \mathbf{B}_m + \mathbf{D}_m.
\]
Finally, it should be noted that the FRF can be computed for a given set of response dofs and a given set of excitation dofs by introducing the appropriate selection matrices \(\mathbf{S}_o\) and \(\mathbf{S}_e\). In this case, the FRF is given by: \[
\mathbf{H}_{oe}(\omega) = \mathbf{S}_o\mathbf{H}(\omega)\mathbf{S}_e.
\]
Note
For the admittance : \(\mathbf{C}_c = [\mathbf{I},\quad \mathbf{0}]\) and \(\mathbf{D}_c = \mathbf{0}\).
For the mobility : \(\mathbf{C}_c = [\mathbf{0}, \quad \mathbf{I}]\) and \(\mathbf{D}_c = \mathbf{0}\).
For the accelerance, \(\mathbf{C}_c = \mathbf{A}_c[m:\text{end}, :]\) and \(\mathbf{D}_c = \mathbf{B}_c[m:\text{end}, :]\) with \(m = N/2\) (this corresponds to a rewrite of the equation of motion).
2.1.1 API
Data types
StateSpaceFRFProblem
StateSpaceFRFProblem(css::ContinuousStateSpace, freq, type, So, Se)
Structure containing the data feeding the direct solver for calculating an FRF
Fields
css: Continuous-time state space model
freq: Frequencies of interest
So: Selection matrix for observation points
Se: Selection matrix for excitation points
Note It is assumed that the output equation is of the form y = So*x
StateSpaceModalFRFProblem
StateSpacemodalFRFProblem(css::ContinuousStateSpace, freq, type, So, Se, n)
Structure containing the data feeding the direct solver for calculating an FRF
Fields
css: Continuous-time state space model
freq: Frequencies of interest
So: Selection matrix for observation points
Se: Selection matrix for excitation points
n: Number of modes to keep in the modal basis
Related function
solve
solve(prob::StateSpaceFRFProblem; type = :dis, ismat = false, progress = true)
solve(prob::StateSpaceModalFRFProblem; type = :dis, ismat = false, progress = true)
Computes the FRF matrix by direct or modal method
Inputs
prob: Structure containing the problem data
type: Type of FRF to compute
:dis: Admittance
:vel: Mobility
:acc: Accelerance
ismat: Return the FRF matrix as a 3D array (default = false)
progress: Show progress bar (default = true)
Output
sol: FRFSolution structure
u: FRF matrix
2.1.2 Example
# Structural parametersM =Diagonal([2., 1.])K = [6. -2.; -2. 4.]ξ =0.05ω, Φ =eigenmode(K, M)C =modal_damping_matrix(M, ω, ξ, Φ)# State-space modelcss =ss_model(K, M, C)# Frequency vectorfreq =0.01:0.001:1.# Problem definition - Case 1 - Directprob_frf =StateSpaceFRFProblem(css, freq)H_direct =solve(prob_frf, ismat =true).u# Problem definition - Case 2 - Modalprob_frf_modal =StateSpaceModalFRFProblem(css, freq)H_modal =solve(prob_frf_modal, ismat =true).u
2.2 Response spectrum
Similarly to the FRF, the response spectrum is a complex function that describes the steady-state response of a system to a sinusoidal input. From the state-space representation, the response spectrum of the system at a given angular frequency \(\omega\) is given by: \[
\mathbf{y}(\omega) = \left[\mathbf{C}_c (j\omega \mathbf{I} - \mathbf{A}_c)^{-1} \mathbf{B}_c + \mathbf{D}_c\right] \mathbf{u}(\omega),
\] or equivalently from the modes of the system: \[
\mathbf{y}(\omega) = \left[\mathbf{C}_m (j\omega \mathbf{I} - \mathbf{\Lambda} )^{-1} \mathbf{B}_m + \mathbf{D}_m\right] \mathbf{u}(\omega).
\]