amfe.structural_dynamics module

Structural Dynamics tools

amfe.structural_dynamics.modal_assurance(U, V)[source]

Compute the Modal Assurance Criterion (MAC) of the vectors stacked in U and V:

\[\begin{split}U = [u_1, \dots, u_n], V = [v_1, \dots, v_n] \\ MAC_{i,j} = \frac{(u_i^Tv_j)^2}{u_i^T u_i \cdot v_j^T v_j}\end{split}\]
Parameters:
  • U (ndarray, shape(N, no_of_modes)) – Array with one set of vectors
  • V (ndarray, shape(N, no_of_modes)) – Array with another set of vectors
Returns:

mac – mac criterion array showing, how the modes coincide. The rows are associated with the vectors in U, the columns with the vectors in V. mac[i,j] gives the squared correlation coefficient of the vecotor U[:,i] and V[:,j].

Return type:

ndarray, shape(no_of_modes, no_of_modes)

References

[R1]Géradin, Michel and Rixen, Daniel: Mechanical Vibrations. John Wiley & Sons, 2014. p.499.
amfe.structural_dynamics.mass_orth(V, M, overwrite=False, niter=2)[source]

Mass-orthogonalize the matrix V with respect to the mass matrix M with a Gram-Schmid-procedure.

Parameters:
  • V (ndarray) – Matrix (e.g. projection basis) containing displacement vectors in the column. Shape is (ndim, no_of_basis_vectors)
  • M (ndarray / sparse matrix.) – Mass matrix. Shape is (ndim, ndim).
  • overwrite (bool) – Flag setting, if matrix V should be overwritten.
  • niter (int) – Number of Gram-Schmid runs for the orthogonalization. As the Gram-Schmid-procedure is not stable, more then one iteration are recommended.
Returns:

V_orth – Mass-orthogonalized basis V

Return type:

ndarray

amfe.structural_dynamics.force_norm(F, K, M, norm='euclidean')[source]

Compute the norm of the given force vector or array

Parameters:
  • F (ndarray, shape (n,) or shape (n, m)) – Array representing the external force in constrained coordinates
  • K (sparse array, shape (n, n)) – Stiffness matrix
  • M (sparse array, shape (n, n)) – Mass matrix
  • norm (str, {'euclidean', 'impedance', 'kinetic'}) –

    norm flag indicating the norm type:

    • ‘euclidean’ : sqrt(F.T @ F)
    • ‘impedance’ : sqrt(F.T @ inv(K) @ F)
    • ‘kinetic’ : sqrt(F.T @ inv(K).T @ M @ inv(K) @ F)
Returns:

norm – norm of the given force vector or array. When F is an array with m columns, norm is a vector with the norm given for every column

Return type:

float or array of shape(m)

amfe.structural_dynamics.rayleigh_coefficients(zeta, omega_1, omega_2)[source]

Compute the coefficients for rayleigh damping such, that the modal damping for the given two eigenfrequencies is zeta.

Parameters:
  • zeta (float) – modal damping for modes 1 and 2
  • omega_1 (float) – first eigenfrequency in [rad/s]
  • omega_2 (float) – second eigenfrequency in [rad/s]
Returns:

  • alpha (float) – rayleigh damping coefficient for the mass matrix
  • beta (float) – rayleigh damping for the stiffness matrix