amfe.reduced_basis module

Reduced basis methods...

amfe.reduced_basis.krylov_subspace(M, K, b, omega=0, no_of_moments=3, mass_orth=True)[source]

Computes the Krylov Subspace associated with the input matrix b at the frequency omega.

Parameters:
  • M (ndarray) – Mass matrix of the system.
  • K (ndarray) – Stiffness matrix of the system.
  • b (ndarray) – input vector of external forcing.
  • omega (float, optional) – frequency for the frequency shift of the stiffness. Default value 0.
  • no_of_moments (int, optional) – number of moments matched. Default value 3.
  • mass_orth (bool, optional) – flag for setting orthogonality of returnd Krylov basis vectors. If True, basis vectors are mass-orthogonal (V.T @ M @ V = eye). If False, basis vectors are orthogonal (V.T @ V = eye)
Returns:

V – Krylov basis where vectors V[:,i] give the basis vectors.

Return type:

ndarray

Examples

TODO

References

amfe.reduced_basis.compute_modes_pardiso(mechanical_system, n=10, u_eq=None, save=False, niter_max=40, rtol=1e-14)[source]

Make a modal analysis using a naive Lanczos iteration.

Parameters:
  • mechanical_system (instance of amfe.MechanicalSystem) – Mechanical system
  • n (int) – number of modes to be computed
  • u_eq (ndarray, optional) – equilibrium position, around which the vibration modes should be computed
  • save (bool, optional) – flat setting, if the modes should be saved in mechanical system. Default value: False
  • niter_max (int, optional) – Maximum number of Lanzcos iterations
Returns:

  • om (ndarray, shape(n)) – eigenfrequencies of the system
  • Phi (ndarray, shape(ndim, n)) – Vibration modes of the system

Note

In comparison to the modal_analysis method, this method uses the Pardiso solver if available and a naive Lanczos iteration. The modal_analysis method uses the arpack solver which is more accurate but takes also much longer for large systems, since the factorization is very inefficient by using superLU.

amfe.reduced_basis.vibration_modes(mechanical_system, n=10, save=False)[source]

Compute the n first vibration modes of the given mechanical system using a power iteration method.

Parameters:
  • mechanical_system (instance of MechanicalSystem) – Mechanical system to be analyzed.
  • n (int) – number of modes to be computed.
  • save (bool) – Flag for saving the modes in mechanical_system for ParaView export. Default: False.
Returns:

  • omega (ndarray) – vector containing the eigenfrequencies of the mechanical system in rad / s.
  • Phi (ndarray) – Array containing the vibration modes. Phi[:,0] is the first vibration mode corresponding to eigenfrequency omega[0]

Examples

Notes

The core command using the ARPACK library is a little bit tricky. One has to use the shift inverted mode for the solution of the mechanical eigenvalue problem with the largest eigenvalues. Generally no convergence is gained when the smallest eigenvalue is to be found.

If the squared eigenvalue omega**2 is negative, as it might happen due to round-off errors with rigid body modes, the negative sign is traveled to the eigenfrequency omega, though this makes physically no sense...

amfe.reduced_basis.craig_bampton(M, K, b, no_of_modes=5, one_basis=True)[source]

Computes the Craig-Bampton basis for the System M and K with the input Matrix b.

Parameters:
  • M (ndarray) – Mass matrix of the system.
  • K (ndarray) – Stiffness matrix of the system.
  • b (ndarray) – Input vector of the system
  • no_of_modes (int, optional) – Number of internal vibration modes for the reduction of the system. Default is 5.
  • one_basis (bool, optional) – Flag for setting, if one Craig-Bampton basis should be returned or if the static and the dynamic basis is chosen separately
Returns:

  • if one_basis=True is chosen
  • V (array) – Basis constisting of static displacement modes and internal vibration modes
  • if one_basis=False is chosen
  • V_static (ndarray) – Static displacement modes corresponding to the input vectors b with V_static[:,i] being the corresponding static displacement vector to b[:,i].
  • V_dynamic (ndarray) – Internal vibration modes with the boundaries fixed.
  • omega (ndarray) – eigenfrequencies of the internal vibration modes.

Examples

TODO

Notes

There is a filter-out command to remove the interface eigenvalues of the system.

References

TODO

amfe.reduced_basis.pod(mechanical_system, n=None)[source]

Compute the POD basis of a mechanical system.

Parameters:
  • mechanical_system (instance of MechanicalSystem) – MechanicalSystem which has run a time simulation and thus displacement fields stored internally.
  • n (int, optional) – Number of POD basis vectors which should be returned. Default is None returning all POD vectors.
Returns:

  • sigma (ndarray) – Array of the singular values.
  • V (ndarray) – Array containing the POD vectors. V[:,0] contains the POD-vector associated with sigma[0] etc.

Examples

TODO

amfe.reduced_basis.modal_derivatives(V, omega, K_func, M, h=1.0, verbose=True, symmetric=True, finite_diff='central')[source]

Compute the basis theta based on real modal derivatives.

Parameters:
  • V (ndarray) – array containing the linear basis
  • omega (ndarray) – eigenfrequencies of the system in rad/s.
  • K_func (function) – function returning the tangential stiffness matrix for a given displacement. Has to work like K = K_func(u).
  • M (ndarray or sparse matrix) – Mass matrix of the system.
  • h (float, optional) – step width for finite difference scheme. Default value is 500 * machine epsilon
  • verbose (bool, optional) – flag for verbosity. Default value: True
  • symmetric (bool, optional) – flag for making the modal derivative matrix theta symmetric. Default is True.
  • finite_diff (str {'central', 'upwind'}) – Method for finite difference scheme. ‘central’ computes the finite difference based on a central difference scheme, ‘upwind’ based on an upwind scheme. Note that the upwind scheme can cause severe distortions of the modal derivative.
Returns:

Theta – three dimensional array of modal derivatives. Theta[:,i,j] contains the modal derivative 1/2 * dx_i / dx_j. The basis Theta is made symmetric, so that Theta[:,i,j] == Theta[:,j,i] if symmetic=True.

Return type:

ndarray

See also

static_correction_theta()
modal derivative with mass neglection.
modal_derivative()
modal derivative for only two vectors.
amfe.reduced_basis.static_derivatives(V, K_func, M=None, omega=0, h=1.0, verbose=True, symmetric=True, finite_diff='central')[source]

Compute the static correction derivatives for the given basis V.

Optionally, a frequency shift can be performed.

Parameters:
  • V (ndarray) – array containing the linear basis
  • K_func (function) – function returning the tangential stiffness matrix for a given displacement. Has to work like K = K_func(u).
  • M (ndarray, optional) – mass matrix. Can be sparse or dense. If None is given, the mass of 0 is assumed. Default value is None.
  • omega (float, optional) – shift frequency. Default value is 0.
  • h (float, optional) – step width for finite difference scheme. Default value is 500 * machine epsilon
  • verbose (bool, optional) – flag for verbosity. Default value: True
  • finite_diff (str {'central', 'forward', backward}) – Method for finite difference scheme. ‘central’ computes the finite difference based on a central difference scheme, ‘forward’ based on an forward scheme etc. Note that the upwind scheme can cause severe distortions of the static correction derivative.
Returns:

Theta – three dimensional array of static corrections derivatives. Theta[:,i,j] contains the static derivative 1/2 * dx_i / dx_j. As the static derivatives are symmetric, Theta[:,i,j] == Theta[:,j,i].

Return type:

ndarray

See also

modal_derivative_theta(), static_correction_derivative()

amfe.reduced_basis.augment_with_derivatives(V, theta, M=None, tol=1e-08, symm=True)[source]

Make a linear basis containing the subspace spanned by V and theta by deflation.

Parameters:
  • V (ndarray) – linear basis
  • theta (ndarray) – third order tensor filled with the modal derivatices associated with V
  • tol (float, optional) – Tolerance for the deflation via SVD. The omitted singular values are at least smaller than the largest one multiplied with tol. Default value: 1E-8.
Returns:

V_ret – linear basis containing the subspace spanned by V and theta.

Return type:

ndarray

amfe.reduced_basis.augment_with_ranked_derivatives(V, Theta, W, n, tol=1e-08, symm=True)[source]

Build a linear basis from a linear basis V and a QM Tensor Theta with a given weighting matrix W.

Parameters:
  • V (ndarray, shape: (ndim, m)) – Linear basis
  • Theta (narray, shape: (ndim, m, m)) – Third order tensor containing the derivatives
  • W (ndarray, shape (m, m)) – Weighting matrix containing the weights
  • n (int) – number of derivative vectors
  • tol (float, optional) – tolerance for the SVD selection. Default value: 1E-8
  • symm (bool, optional) – flag if weighting matrix and Theta is symmetric. If true, only the upper half of Theta and W is evaluated
Returns:

V_red

Return type:

ndarray, shape: (ndim, n + m)

amfe.reduced_basis.ranking_of_weighting_matrix(W, symm=True)[source]

Return the indices of a weighting matrix W in decreasing order

Parameters:
  • W (ndarray, shape: (n, m)) – weighting matrix
  • symm (bool, optional) – flag indicating symmetry
Returns:

idxs – array of indices, where i, j = idxs[k] are the indices of the k-largest value in W

Return type:

ndarray, shape (n*m, 2)