Description
Modal analysis and substructuring.
Using this module, you can perform modal analysis directly in Chrono. This can be useful for computing natural frequencies, for computing stability (complex eigenvalue analysis, with damping matrices), or for substructuring (where subassemblies are replaced by modal bodies).
For additional information, see:
- the installation guide
- the tutorials
Modules | |
Run-time visualization | |
Namespaces | |
chrono::modal | |
Namespace with classes for the modal module. | |
Classes | |
class | chrono::modal::ChGeneralizedEigenvalueSolver< ScalarT > |
Base interface class for generalized eigenvalue solvers A*x = lambda*B*x. More... | |
class | chrono::modal::callback_Ax |
Generic A*x callback. More... | |
class | chrono::modal::callback_Ax_sparse_shiftinvert |
The callback to be used for "A*x" where for shift&invert is: A = (As - sigma Bs)/Bs , so A*x = (As - sigma Bs)/(Bs*x), just like a linear system with coefficient matrix (As - sigma Bs) and known rhs Bs*x. More... | |
class | chrono::modal::callback_Ax_sparse_complexshiftinvert |
The callback to be used for "A*x" where for shift&invert is: A = (As - sigma Bs)/Bs , with COMPLEX sigma shift, so A*x = (As - sigma Bs)/(Bs*x), just like a linear system with coefficient matrix (As - sigma Bs) and known rhs Bs*x. More... | |
class | chrono::modal::ChKrylovSchurEig |
Compute (complex) eigenvalues and eigenvectors using the Krylov-Schur algorithm. More... | |
class | chrono::modal::ChModalAssembly |
Class for assemblies of items, for example ChBody, ChLink, ChMesh, etc. More... | |
class | chrono::modal::ChModalDamping |
Base class for damping models of modal reduced assemblies. More... | |
class | chrono::modal::ChModalDampingNone |
Class for no damping model. More... | |
class | chrono::modal::ChModalDampingRayleigh |
Class for simple Rayleigh damping model R^ = alpha*M^ + beta*K^ where M^ and K^ are the reduced matrices, both for boundary nodes and modal coords. More... | |
class | chrono::modal::ChModalDampingFactorRmm |
Class for setting the damping via N damping factors z_i of the internal mode coordinates. More... | |
class | chrono::modal::ChModalDampingFactorRayleigh |
Class for setting the damping via N damping factors z_i of the internal mode coordinates and alpha-beta Rayleigh damping for the boundary nodes, assuming R^ = [Rbb Rbm ] [Rmb Rmm ] with Rmm=diag { 2 z_1 w_1, 2 z_2 w_2, ..., 2 z_i w_i }, Rbb= alpha*Mbb + beta*Kbb, Rbm = 0, Rmb = 0. More... | |
class | chrono::modal::ChModalDampingFactorAssembly |
Class for setting the damping via N damping factors z_i for all the modes of the subassembly, where assembly n.modes = (boundary coords+internal modes) R^ = V'^-1 * Dd * V^-1 with Dd=diag { 2 z_1 w_1, 2 z_2 w_2, ..., 2 z_i w_i }, and V = eigenvectors of (M^, K^). More... | |
class | chrono::modal::ChModalDampingReductionR |
Class for damping as reduction of the original damping matrix via the eigenvectors of the undamped assembly, i.e. More... | |
class | chrono::modal::ChModalDampingCustom |
Class for damping defined with an user-defined matrix that could be obtained via external tools such as Matlab or FEA. More... | |
class | chrono::modal::ChModalSolver |
Base class for modal solvers. More... | |
class | chrono::modal::ChModalSolverDamped |
Modal solver for damped systems of the form (-w^2*M + i*w*R + K)*x = 0 s.t. More... | |
class | chrono::modal::ChModalSolverUndamped< EigenvalueSolverType > |
Modal solver for undamped systems (-w^2*M+K)*x = 0 s.t. More... | |
class | chrono::modal::ChSymGenEigenvalueSolver |
Base interface class for iterative eigenvalue solvers for generalized problem with real symmetric matrices. More... | |
class | chrono::modal::ChSymGenEigenvalueSolverKrylovSchur |
Generalized iterative eigenvalue solver implementing Krylov-Schur shift-and-invert method for real symmetric matrices. More... | |
class | chrono::modal::ChSymGenEigenvalueSolverLanczos |
Generalized iterative eigenvalue solver implementing Lanczos shift-and-invert method for real symmetric matrices. More... | |
class | chrono::modal::ChUnsymGenEigenvalueSolver |
Base interface class for iterative eigenvalue solvers for generalized problem with real generic matrices. More... | |
class | chrono::modal::ChUnsymGenEigenvalueSolverKrylovSchur |
Generalized iterative eigenvalue solver implementing Krylov-Schur shift-and-invert method for real generic matrices. More... | |
Functions | |
template<typename EigSolverType > | |
int | chrono::modal::Solve (EigSolverType &eig_solver, ChSparseMatrix &A, ChSparseMatrix &B, ChMatrixDynamic< typename EigSolverType::ScalarType > &eigvects, ChVectorDynamic< typename EigSolverType::ScalarType > &eigvals, const std::list< std::pair< int, typename EigSolverType::ScalarType >> &eig_requests, bool uniquify=true, int eigvects_clipping_length=0) |
Helper function to solve any kind of generalized eigenvalue problem even with multiple requests. More... | |
void | chrono::modal::BuildUndampedEigenProblemMatrices (ChAssembly &assembly, ChSystemDescriptor &temp_descriptor, ChSparseMatrix &A, ChSparseMatrix &B, int n_vars) |
Partially build the undamped eigenvalue problem matrices A and B from a given ChAssembly. More... | |
template<typename EigensolverType > | |
void | chrono::modal::ChModalAssembly::DoModalReduction (const ChModalSolverUndamped< EigensolverType > &modal_solver, const ChModalDamping &damping_model=ChModalDampingNone()) |
Perform modal reduction on this modal assembly, from the current "full" ("boundary"+"internal") assembly. More... | |
template<typename EigensolverType > | |
void | chrono::modal::ChModalAssembly::DoModalReduction (ChSparseMatrix &full_M, ChSparseMatrix &full_K, ChSparseMatrix &full_Cq, const ChModalSolverUndamped< EigensolverType > &modal_solver, const ChModalDamping &damping_model=ChModalDampingNone()) |
Perform modal reduction on this modal assembly that contains only the "boundary" nodes, whereas the "internal" nodes have been modeled only in an external FEA software with the full ("boundary"+"internal") modes. More... | |
int | chrono::modal::ChModalSolverUndamped< EigenvalueSolverType >::Solve (const ChAssembly &assembly, ChMatrixDynamic< ScalarType > &eigvects, ChVectorDynamic< ScalarType > &eigvals, ChVectorDynamic< double > &freq) const |
Solve the constrained eigenvalue problem retrieving it from the ChAssembly. More... | |
int | chrono::modal::ChModalSolverUndamped< EigenvalueSolverType >::Solve (const ChSparseMatrix &K, const ChSparseMatrix &M, const ChSparseMatrix &Cq, ChMatrixDynamic< ScalarType > &eigvects, ChVectorDynamic< ScalarType > &eigvals, ChVectorDynamic< double > &freq) const |
Solve the constrained eigenvalue problem setting it up from individual matrices. More... | |
Function Documentation
◆ BuildUndampedEigenProblemMatrices()
void ChApiModal chrono::modal::BuildUndampedEigenProblemMatrices | ( | ChAssembly & | assembly, |
ChSystemDescriptor & | temp_descriptor, | ||
ChSparseMatrix & | A, | ||
ChSparseMatrix & | B, | ||
int | n_vars | ||
) |
Partially build the undamped eigenvalue problem matrices A and B from a given ChAssembly.
Utility function to build the A and B matrices for the undamped eigenvalue problem.
WARNING: Cq and Cq' signs are not flipped here: the user is expected to flip it during scaling (if any). The finale shape should be: A = [ -K -Cq' ] [ -Cq 0 ]
B = [ M 0 ] [ 0 0 ]
◆ DoModalReduction() [1/2]
void chrono::modal::ChModalAssembly::DoModalReduction | ( | ChSparseMatrix & | full_M, |
ChSparseMatrix & | full_K, | ||
ChSparseMatrix & | full_Cq, | ||
const ChModalSolverUndamped< EigensolverType > & | modal_solver, | ||
const ChModalDamping & | damping_model = ChModalDampingNone() |
||
) |
Perform modal reduction on this modal assembly that contains only the "boundary" nodes, whereas the "internal" nodes have been modeled only in an external FEA software with the full ("boundary"+"internal") modes.
- with an external FEA software, the full assembly is modeled with "boundary"+"internal" nodes.
- with an external FEA software, the M mass matrix and the K stiffness matrix are saved to disk.
- in Chrono, M and K and Cq constraint jacobians (if any) are load from disk and stored in ChSparseMatrix objects.
- in Chrono, only boundary nodes are added to a ChModalAssembly.
- in Chrono, run this function passing such M and K matrices: a modal analysis will be done on K and M Note that the size of M (and K) must be at least > m_num_coords_vel_boundary.
- Parameters
-
full_M mass matrix of the full assembly (boundary+internal) full_K stiffness matrix of the full assembly (boundary+internal) full_Cq constraint jacobian matrix of the full assembly (boundary+internal) modal_solver settings for the modal analysis, such as the number of modes to extract damping_model damping model
◆ DoModalReduction() [2/2]
void chrono::modal::ChModalAssembly::DoModalReduction | ( | const ChModalSolverUndamped< EigensolverType > & | modal_solver, |
const ChModalDamping & | damping_model = ChModalDampingNone() |
||
) |
Perform modal reduction on this modal assembly, from the current "full" ("boundary"+"internal") assembly.
- An undamped modal analysis will be done on the full assembly, followed by a modal reduction transformation.
- The "boundary" nodes will be retained.
- The "internal" nodes will be replaced by n_modes modal coordinates.
◆ Solve() [1/3]
int chrono::modal::ChModalSolverUndamped< EigenvalueSolverType >::Solve | ( | const ChAssembly & | assembly, |
ChMatrixDynamic< ScalarType > & | eigvects, | ||
ChVectorDynamic< ScalarType > & | eigvals, | ||
ChVectorDynamic< double > & | freq | ||
) | const |
Solve the constrained eigenvalue problem retrieving it from the ChAssembly.
Only the position part of the eigenvectors is returned, unless SetClipPositionCoords(false) is called.
◆ Solve() [2/3]
int chrono::modal::ChModalSolverUndamped< EigenvalueSolverType >::Solve | ( | const ChSparseMatrix & | K, |
const ChSparseMatrix & | M, | ||
const ChSparseMatrix & | Cq, | ||
ChMatrixDynamic< ScalarType > & | eigvects, | ||
ChVectorDynamic< ScalarType > & | eigvals, | ||
ChVectorDynamic< double > & | freq | ||
) | const |
Solve the constrained eigenvalue problem setting it up from individual matrices.
Only the position part of the eigenvectors is returned, unless SetClipPositionCoords(false) is called.
◆ Solve() [3/3]
int chrono::modal::Solve | ( | EigSolverType & | eig_solver, |
ChSparseMatrix & | A, | ||
ChSparseMatrix & | B, | ||
ChMatrixDynamic< typename EigSolverType::ScalarType > & | eigvects, | ||
ChVectorDynamic< typename EigSolverType::ScalarType > & | eigvals, | ||
const std::list< std::pair< int, typename EigSolverType::ScalarType >> & | eig_requests, | ||
bool | uniquify = true , |
||
int | eigvects_clipping_length = 0 |
||
) |
Helper function to solve any kind of generalized eigenvalue problem even with multiple requests.
The eigenvectors are also filtered to avoid duplicates with the same frequency. This means that only one of the complex eigenvectors that come in conjugate pairs is stored. The eigrequest argument is a list of pairs, where the first element is the number of modes to be found, and the second element is the shift to apply for that specific search.