Multiple Degree of Freedom Systems (vibration_toolbox.mdof
)¶
mdof
– Multiple Degree of Freedom Functions¶
Multiple Degree of Freedom Analysis Tools.
-
mdof.
modes_system
(M, K, C=None)[source]¶ Natural frequencies, damping ratios, and mode shapes of MDOF system. This function will return the natural frequencies (wn), the damped natural frequencies (wd), the damping ratios (zeta), the right eigenvectors (X) and the left eigenvectors (Y) for a system defined by M, K and C. If the dampind matrix ‘C’ is none or if the damping is proportional, wd and zeta will be none and X and Y will be equal.
- Parameters
- M: array
Mass matrix
- K: array
Stiffness matrix
- C: array
Damping matrix
- Returns
- wn: array
The natural frequencies of the system
- wd: array
The damped natural frequencies of the system
- zeta: array
The damping ratios
- X: array
The right eigenvectors
- Y: array
The left eigenvectors
Examples
>>> # This test has been moved to tests_mdof >>> M = np.array([[1, 0], ... [0, 1]]) >>> K = np.array([[2, -1], ... [-1, 6]]) >>> C = np.array([[0.3, -0.02], ... [-0.02, 0.1]]) >>> wn, wd, zeta, X, Y = modes_system(M, K, C) # doctest: +SKIP Damping is non-proportional, eigenvectors are complex. >>> wn # doctest: +SKIP array([1.33, 2.5 , 1.33, 2.5 ]) >>> wd # doctest: +SKIP array([1.32, 2.5 , 1.32, 2.5 ]) >>> zeta # doctest: +SKIP array([0.11, 0.02, 0.11, 0.02]) >>> X # doctest: +SKIP array([[-0.06-0.58j, -0.01+0.08j, -0.06+0.58j, -0.01-0.08j], [-0. -0.14j, -0.01-0.36j, -0. +0.14j, -0.01+0.36j], [ 0.78+0.j , -0.21-0.03j, 0.78-0.j , -0.21+0.03j], [ 0.18+0.01j, 0.9 +0.j , 0.18-0.01j, 0.9 -0.j ]]) >>> Y # doctest: +SKIP array([[ 0.02+0.82j, 0.01-0.31j, 0.02-0.82j, 0.01+0.31j], [-0.05+0.18j, 0.01+1.31j, -0.05-0.18j, 0.01-1.31j], [ 0.61+0.06j, -0.12-0.02j, 0.61-0.06j, -0.12+0.02j], [ 0.14+0.03j, 0.53+0.j , 0.14-0.03j, 0.53-0.j ]]) >>> C = 0.2*K # with proportional damping >>> wn, wd, zeta, X, Y = modes_system(M, K, C) # doctest: +SKIP Damping is proportional or zero, eigenvectors are real >>> X # doctest: +SKIP array([[-0.97, 0.23], [-0.23, -0.97]])
-
mdof.
modes_system_undamped
(M, K)[source]¶ Return eigensolution of multiple DOF system.
Returns the natural frequencies (w), eigenvectors (P), mode shapes (S) and the modal transformation matrix S for an undamped system.
See Notes for explanation of the underlying math.
- Parameters
- M: float array
Mass matrix
- K: float array
Stiffness matrix
- Returns
- w: float array
The natural frequencies of the system
- P: float array
The eigenvectors of the system.
- S: float array
The mass-normalized mode shapes of the system.
- Sinv: float array
The modal transformation matrix S^-1(takes x -> r(modal coordinates))
Notes
Given \(M\ddot{x}(t)+Kx(t)=0\), with mode shapes \(u\), the matrix of mode shapes \(S=[u_1 u_1 \ldots]\) can be created. If the modal coordinates are the vector \(r(t)\). The modal transformation separates space and time from \(x(t)\) such that \(x(t)=S r(t)\). Substituting into the governing equation:
\(MS\ddot{r}(t)+KSr(t)=0\)
Premultiplying by \(S^T\)
\(S^TMS\ddot{r}(t)+S^TKSr(t)=0\)
The matrices \(S^TMS\) and \(S^TKS\) will be diagonalized by this process (\(u_i\) are the eigenvectors of \(M^{-1}K\)).
If scaled properly (mass normalized so \(u_i^TMu_i=1\)) then \(S^TMS=I\) and \(S^TKS=\Omega^2\) where \(\Omega^2\) is a diagonal matrix of the natural frequencies squared in radians per second.
Further, inverses are unstable so the better way to solve linear equations is with Gauss elimination.
\(AB=C\) given known \(A\) and \(C\) is solved using la.solve(A, C, assume_a=’pos’).
\(BA=C\) given known \(A\) and \(C\) is solved by first transposing the equation to \(A^TB^T=C^T\), then solving for \(C^T\). The resulting command is la.solve(A.T, C.T, assume_a=’pos’).T
Examples
>>> M = np.array([[4, 0, 0], ... [0, 4, 0], ... [0, 0, 4]]) >>> K = np.array([[8, -4, 0], ... [-4, 8, -4], ... [0, -4, 4]]) >>> w, P, S, Sinv = modes_system_undamped(M, K) >>> w # doctest: +SKIP array([0.45, 1.25, 1.8 ]) >>> S array([[ 0.16, -0.37, -0.3 ], [ 0.3 , -0.16, 0.37], [ 0.37, 0.3 , -0.16]])
-
mdof.
response_system
(M, C, K, F, x0, v0, t)[source]¶ Returns system response given the initial displacement vector ‘X0’, initial velocity vector ‘V0’, the mass matrix ‘M’, the stiffness matrix ‘M’, and the damping matrix ‘C’ and force ‘F’. T is a row vector of evenly spaced times. F is a matrix of forces over time, each column corresponding to the corresponding column of T, each row corresponding to the same numbered DOF.
- Parameters
- M: array
Mass matrix
- K: array
Stiffness matrix
- C: array
Damping matrix
- x0: array
Array with displacement initial conditions
- v0: array
Array with velocity initial conditions
- t: array
Array withe evenly spaced times
- Returns
- Tarray
Time values for the output.
- youtarray
System response.
- xoutarray
Time evolution of the state vector.
Examples
>>> M = np.array([[9, 0], ... [0, 1]]) >>> K = np.array([[27, -3], ... [-3, 3]]) >>> C = K/10 >>> x0 = np.array([0, 1]) >>> v0 = np.array([1, 0]) >>> t = np.linspace(0, 10, 100) >>> F = np.vstack([0*t, ... 3*np.cos(2*t)]) >>> tou, yout, xout = response_system(M, C, K, F, x0, v0, t) >>> tou[:10] # doctest: +SKIP array([0. , 0.1 , 0.2 , 0.3 , 0.4 , 0.51, 0.61, 0.71, 0.81, 0.91])
>>> yout[:10] array([[ 0. , 1. , 1. , 0. ], [ 0.1 , 1. , 0.99, 0.04], [ 0.2 , 1.01, 0.95, 0.1 ], [ 0.29, 1.02, 0.88, 0.15], [ 0.38, 1.04, 0.79, 0.19], [ 0.45, 1.06, 0.68, 0.2 ], [ 0.51, 1.08, 0.55, 0.17], [ 0.56, 1.09, 0.41, 0.09], [ 0.59, 1.09, 0.26, -0.04], [ 0.61, 1.08, 0.11, -0.22]])
>>> xout[:10] array([[ 0. , 1. , 1. , 0. ], [ 0.1 , 1. , 0.99, 0.04], [ 0.2 , 1.01, 0.95, 0.1 ], [ 0.29, 1.02, 0.88, 0.15], [ 0.38, 1.04, 0.79, 0.19], [ 0.45, 1.06, 0.68, 0.2 ], [ 0.51, 1.08, 0.55, 0.17], [ 0.56, 1.09, 0.41, 0.09], [ 0.59, 1.09, 0.26, -0.04], [ 0.61, 1.08, 0.11, -0.22]])
-
mdof.
response_system_undamped
(M, K, x0, v0, max_time)[source]¶ This function calculates the time response for an undamped system and returns the vector (state-space) X. The n first rows contain the displacement (x) and the n last rows contain velocity (v) for each coordinate. Each column is related to a time-step. The time array is also returned.
- Parameters
- M: array
Mass matrix
- K: array
Stiffness matrix
- x0: array
Array with displacement initial conditions
- v0: array
Array with velocity initial conditions
- max_time: float
End time
- Returns
- t: array
Array with the time
- X: array
The state-space vector for each time
Examples
>>> M = np.array([[1, 0], ... [0, 4]]) >>> K = np.array([[12, -2], ... [-2, 12]]) >>> x0 = np.array([1, 1]) >>> v0 = np.array([0, 0]) >>> max_time = 10 >>> t, X = response_system_undamped(M, K, x0, v0, max_time) >>> # first column is the initial conditions [x1, x2, v1, v2] >>> X[:, 0] # doctest: +SKIP array([1., 1., 0., 0.]) >>> X[:, 1] # displacement and velocities after delta t array([ 1. , 1. , -0.04, -0.01])