Modal¶
- class sharpy.solvers.modal.Modal[source]¶
Modal
solver class, inherited fromBaseSolver
Extracts the
M
,K
andC
matrices from theFortran
library for the beam. Depending on the choice of modal projection, these may or may not be transformed to a state-space form to compute the eigenvalues and mode shapes of the structure.The settings that this solver accepts are given by a dictionary, with the following key-value pairs:
Name
Type
Description
Default
print_info
bool
Write status to screen
True
rigid_body_modes
bool
Write modes with rigid body mode shapes
False
use_undamped_modes
bool
Project the modes onto undamped mode shapes
True
NumLambda
int
Number of modes to retain
20
write_modes_vtk
bool
Write Paraview files with mode shapes
True
print_matrices
bool
Write M, C and K matrices to file
False
save_data
bool
Write mode shapes, frequencies and damping to file
True
continuous_eigenvalues
bool
Use continuous time eigenvalues
False
dt
float
Time step to compute discrete time eigenvalues
0
delta_curved
float
Threshold for linear expressions in rotation formulas
0.01
plot_eigenvalues
bool
Plot to screen root locus diagram
False
max_rotation_deg
float
Scale mode shape to have specified maximum rotation
15.0
max_displacement
float
Scale mode shape to have specified maximum displacement
0.15
use_custom_timestep
int
If > -1, it will use that time step geometry for calculating the modes
-1
rigid_modes_ppal_axes
bool
Modify the ridid body modes such that they are defined wrt to the CG and aligned with the principal axes of inertia
False
rigid_modes_cg
bool
Not implemente yet
False
- free_free_modes(phi, M)[source]¶
Warning
This function is deprecated. See
free_modes_principal_axes()
for a transformation to the CG and with respect to the principal axes of inertia.Returns the rigid body modes defined with respect to the centre of gravity
The transformation from the modes defined at the FoR A origin, \(\boldsymbol{\Phi}\), to the modes defined using the centre of gravity as a reference is
\[\boldsymbol{\Phi}_{rr,CG}|_{TRA} = \boldsymbol{\Phi}_{RR}|_{TRA} + \tilde{\mathbf{r}}_{CG} \boldsymbol{\Phi}_{RR}|_{ROT}\]\[\boldsymbol{\Phi}_{rr,CG}|_{ROT} = \boldsymbol{\Phi}_{RR}|_{ROT}\]- Returns
Transformed eigenvectors
- Return type
(np.array)
- run(**kwargs)[source]¶
Extracts the eigenvalues and eigenvectors of the clamped structure.
If
use_undamped_modes == True
then the free vibration modes of the clamped structure are found solving:\[\mathbf{M\,\ddot{\eta}} + \mathbf{K\,\eta} = 0\]that flows down to solving the non-trivial solutions to:
\[(-\omega_n^2\,\mathbf{M} + \mathbf{K})\mathbf{\Phi} = 0\]On the other hand, if the damped modes are chosen because the system has damping, the free vibration modes are found solving the equation of motion of the form:
\[\mathbf{M\,\ddot{\eta}} + \mathbf{C\,\dot{\eta}} + \mathbf{K\,\eta} = 0\]which can be written in state space form, with the state vector \(\mathbf{x} = [\eta^T,\,\dot{\eta}^T]^T\) as
\[\begin{split}\mathbf{\dot{x}} = \begin{bmatrix} 0 & \mathbf{I} \\ -\mathbf{M^{-1}K} & -\mathbf{M^{-1}C} \end{bmatrix} \mathbf{x}\end{split}\]and therefore the mode shapes and frequencies correspond to the solution of the eigenvalue problem
\[\mathbf{A\,\Phi} = \mathbf{\Lambda\,\Phi}.\]From the eigenvalues, the following system characteristics are provided:
Natural Frequency: \(\omega_n = |\lambda|\)
Damped natural frequency: \(\omega_d = \text{Im}(\lambda) = \omega_n \sqrt{1-\zeta^2}\)
Damping ratio: \(\zeta = -\frac{\text{Re}(\lambda)}{\omega_n}\)
In addition to the above, the modal output dictionary includes the following:
M
: Tangent mass matrixC
: Tangent damping matrixK
: Tangent stiffness matrixCcut
: Modal damping matrix \(\mathbf{C}_m = \mathbf{\Phi}^T\mathbf{C}\mathbf{\Phi}\)Kin_damp
: Forces gain matrix (when damped): \(K_{in} = \mathbf{\Phi}_L^T \mathbf{M}^{-1}\)eigenvectors
: Right eigenvectorseigenvectors_left
: Left eigenvectors given when the system is damped
- Returns
updated data object with modal analysis as part of the last structural time step.
- Return type
PreSharpy