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: |
|
---|---|
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: |
|
---|---|
Returns: |
|
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: |
|
---|---|
Returns: |
|
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: |
|
---|---|
Returns: |
|
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: |
|
---|---|
Returns: |
|
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: |
|
---|---|
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()
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: |
|
---|---|
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: |
|
---|---|
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: |
|
---|---|
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: |
|
---|---|
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) |