From a general perspective, the frequency-domain solvers are used to analyze the dynamic response of structures subjected to harmonic loads. In the frequency domain, Mdof systems are described at a given angular frequency \(\omega\) by the following equation: \[
\left(\mathbf{K} + j\omega\mathbf{C} - \omega^2\mathbf{M}\right) \mathbf{x}(\omega) = \mathbf{F}(\omega).
\]
1 Frequency response function
The Frequency response function (FRF) is computed by assuming a unit harmonic load at the \(i\)-th degree of freedom (DOF) and zero loads at all other DOFs, meaning that \(\mathbf{F}(\omega) = \mathbf{I}\).
Consequently, broadly speaking, the FRF is defined by: \[
\mathbf{H}(\omega) = \left(\mathbf{K} + j\omega\mathbf{C} - \omega^2\mathbf{M}\right)^{-1}.
\] This method is known as the direct method.
Alternatively, the FRF can be computed using the modal method, which is based on the eigenvalue decomposition of the system matrices. The modal method is more efficient for large systems, as it reduces the problem to a smaller number of modes. In this case and assuming a modal damping ratio, the FRF is computed as: \[
\mathbf{H}(\omega) = \mathbf{\Phi}\left(\mathbf{\Omega}^2 - \omega^2\mathbf{I} + 2j\omega\mathbf{\Xi}\right)^{-1}\mathbf{\Phi}^T,
\] where \(\mathbf{\Phi}\) is the mass-normalized mode shapes matrix and \(\mathbf{\Omega}^2 = \text{diag}(\omega_1^2, \dots, \omega_n^2)\) and \(\mathbf{\Xi} = \text{diag}(\xi_1\omega_1, \dots, \xi_n\omega_n)\).
Note
As explained in State-Space solvers - Section 2.1, 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.
\]
1.1 API
Data types
DirectFRFProblem
DirectFRFProblem(K, M, C, freq, So, Se)
Structure containing the data feeding the direct solver for calculating an FRF
Fields
K::AbstractMatrix: Stiffness matrix
M::AbstractMatrix: Mass matrix
C::AbstractMatrix: Damping matrix
freq::AbtractRange: Frequencies of interest
So::AbstractMatrix: Selection matrix for observation points
Se::AbstractMatrix: Selection matrix for excitation points
ModalFRFProblem
ModalFRFProblem(ωn, ξn, freq, ϕo, ϕe)
Structure containing the data feeding the modal solver for calculating an FRF
Fields
ωn::Vector{Real}: Natural angular frequencies
ξn::Vector{Real}: Modal damping ratios
freq::AbstractRange: Frequencies of interest
ϕe::AbstractMatrix: Mode shapes at excitation points
ϕo::AbstractMatrix: Mode shapes at observation points
Note The mode shapes must be mass-normalized
Related function
solve
solve(prob::DirectFRF; type = :dis, ismat = false, progress = true)
solve(prob::ModalFRFProblem; type = :dis, ismat = false, progress = true)
Computes the FRF matrix by direct or modal approach
Inputs
prob: Structure containing the problem data
type: Type of FRF to compute
:dis: Admittance (default)
:vel: Mobility
:acc: Accelerance
ismat: Return the FRF matrix as a 3D array (default = false)
Similarly to the FRF, the response spectrum is a complex function that describes the steady-state response of a system to a sinusoidal input. At a given angular frequency \(\omega\), it is expressed as : \[
\mathbf{x}(\omega) = \left(\mathbf{K} + j\omega\mathbf{C} - \omega^2\mathbf{M}\right)^{-1} \mathbf{F}(\omega).
\] or equivalently from the modes of the system: \[
\mathbf{x}(\omega) = \mathbf{\Phi}\left(\mathbf{\Omega}^2 - \omega^2\mathbf{I} + 2j\omega\mathbf{\Omega}\right)^{-1}\mathbf{L},
\] where \(\mathbf{L} = \mathbf{\Phi}^\mathsf{T}\mathbf{F}(\omega)\) is the modal participation factor matrix1.
2.1 API
Data types
DirectFreqProblem
DirectFreqProblem(K, M, C, F, freq, So)
Structure containing the data feeding the direct solver for calculating the modal frequencies
Fields
K::AbstractMatrix: Stiffness matrix
M::AbstractMatrix: Mass matrix
C::AbstractMatrix: Damping matrix
F::AbstractMatrix: Force matrix
freq::AbstractRange: Frequencies of interest
So::AbstractMatrix: Selection matrix for observation points
ModalFreqProblem
ModalFreqProblem(ωn, ξn, Fn, freq, ϕo)
Structure containing the data feeding the modal solver for calculating the frequency response by modal approach
Fields
ωn::Vector{Real}: Natural angular frequencies
ξn::Vector{Real}: Modal damping ratios
Ln::AbstractMatrix: Modal participation factors
freq::AbstractRange: Frequencies of interest
ϕo::Matrix{Real}: Mode shapes at observation points
Related function
solve
solve(prob::DirectFreqProblem, type = :dis, progress = true)
solve(prob::ModalFreqProblem, type = :dis, progress = true)
Computes the frequency response by direct or modal approach
Inputs
prob: Structure containing the problem data
type: Type of response to compute
:dis: Displacement (default)
:vel: Velocity
:acc: Acceleration
progress: Show progress bar (default = true)
Output
sol: Solution of the problem
u: Response spectrum matrix
2.2 Example
# Structural parametersM =Diagonal([2., 1.])K = [6. -2.; -2. 4.]ξ =0.05ω, Φ =eigenmode(K, M)C =modal_damping_matrix(M, ω, ξ, Φ)# Frequency vectorfreq =0.01:0.001:1.# Force matrixF =zeros(2, length(freq))F[1, :] .=10.# Problem definition - Case 1 - Directprob_freq =DirectFreqProblem(K, M, C, F, freq)y_freq =solve(prob_freq).u# Problem definition - Case 2 - Modalprob_freq_modal =ModalFreqProblem(ω, ξ, Φ'*F, freq, Φ)y_freq_modal =solve(prob_freq_modal).u
Footnotes
The expression of the modal participation factor matrix is the same as the modal force matrix, because the mode shapes are assumed to be mass-normalized.↩︎