chrono::modal::ChModalSolverDamped Class Reference

Description

Modal solver for damped systems of the form (-w^2*M + i*w*R + K)*x = 0 s.t.

Cq*x = 0. with complex w (where w.length() = undamped nat.freq)

#include <ChModalSolverDamped.h>

Inheritance diagram for chrono::modal::ChModalSolverDamped:
Collaboration diagram for chrono::modal::ChModalSolverDamped:

Public Member Functions

 ChModalSolverDamped (int n_lower_modes, double base_freq, bool scaleCq=true, bool verbose=false, std::shared_ptr< ChUnsymGenEigenvalueSolver > solver=chrono_types::make_shared< ChUnsymGenEigenvalueSolverKrylovSchur >())
 Creates a modal solver for the damped case. More...
 
 ChModalSolverDamped (std::vector< ChFreqSpan > freq_spans, bool scaleCq=true, bool verbose=false, std::shared_ptr< ChUnsymGenEigenvalueSolver > solver=chrono_types::make_shared< ChUnsymGenEigenvalueSolverKrylovSchur >())
 Creates a modal solver for the damped case. More...
 
virtual int Solve (const ChAssembly &assembly, ChMatrixDynamic< std::complex< double >> &eigvects, ChVectorDynamic< std::complex< double >> &eigvals, ChVectorDynamic< double > &freq, ChVectorDynamic< double > &damp_ratios) const
 Solve the constrained eigenvalue problem retrieving it from the ChAssembly. More...
 
virtual int Solve (const ChSparseMatrix &K, const ChSparseMatrix &R, const ChSparseMatrix &M, const ChSparseMatrix &Cq, ChMatrixDynamic< std::complex< double >> &eigvects, ChVectorDynamic< std::complex< double >> &eigvals, ChVectorDynamic< double > &freq, ChVectorDynamic< double > &damp_ratios) const
 Solve the constrained eigenvalue problem setting it up from individual matrices. More...
 
std::shared_ptr< ChUnsymGenEigenvalueSolverGetEigenSolver () const
 Get the inner eigensolver.
 
- Public Member Functions inherited from chrono::modal::ChModalSolver
 ChModalSolver (int n_lower_modes, double base_freq, bool scaleCq, bool verbose)
 Creates a modal solver. More...
 
 ChModalSolver (std::vector< ChFreqSpan > freq_spans, bool scaleCq, bool verbose)
 Creates a modal solver. More...
 
int GetNumRequestedModes () const
 Get the total number of requested modes.
 
void SetClipPositionCoords (bool val)
 Clip the eigenvectors to only the position coordinates.
 
const std::vector< ChFreqSpan > & GetFrequencySpans () const
 Get the set of frequency spans for which modes are requested.
 
double GetTimeMatrixAssembly () const
 Get cumulative time for matrix assembly.
 
double GetTimeEigenSolver () const
 Get cumulative time eigensolver solution.
 
double GetTimeSolutionPostProcessing () const
 Get cumulative time for post-solver solution postprocessing.
 

Protected Attributes

std::shared_ptr< ChUnsymGenEigenvalueSolverm_solver
 
- Protected Attributes inherited from chrono::modal::ChModalSolver
ChTimer m_timer_matrix_assembly
 timer for matrix assembly
 
ChTimer m_timer_eigen_solver
 timer for eigensolver solution
 
ChTimer m_timer_solution_postprocessing
 timer for conversion of eigensolver solution
 
std::vector< ChFreqSpan > m_freq_spans
 frequency spans for which modes are requested
 
bool m_clip_position_coords
 store only the part of each eigenvector that refers to the position coordinates More...
 
bool m_scaleCq = true
 if true, the Cq matrix is scaled to improve conditioning
 
bool m_verbose = false
 if true, additional information is printed during the solution process
 

Constructor & Destructor Documentation

◆ ChModalSolverDamped() [1/2]

chrono::modal::ChModalSolverDamped::ChModalSolverDamped ( int  n_lower_modes,
double  base_freq,
bool  scaleCq = true,
bool  verbose = false,
std::shared_ptr< ChUnsymGenEigenvalueSolver solver = chrono_types::make_shared<ChUnsymGenEigenvalueSolverKrylovSchur>() 
)
inline

Creates a modal solver for the damped case.

n_lower_modes number of lowest modes to be found. base_freq frequency around which the modes will be found; higher values can help finding eigenvalues for ill-conditioned problems. scaleCq if true, the Cq matrix is scaled to improve conditioning. verbose if true, additional information is printed during the solution process. solver the inner eigensolver to be used.

◆ ChModalSolverDamped() [2/2]

chrono::modal::ChModalSolverDamped::ChModalSolverDamped ( std::vector< ChFreqSpan >  freq_spans,
bool  scaleCq = true,
bool  verbose = false,
std::shared_ptr< ChUnsymGenEigenvalueSolver solver = chrono_types::make_shared<ChUnsymGenEigenvalueSolverKrylovSchur>() 
)
inline

Creates a modal solver for the damped case.

freq_spans pair of number of modes and frequency around which the modes will be found. ill-conditioned problems. scaleCq if true, the Cq matrix is scaled to improve conditioning. verbose if true, additional information is printed during the solution process. solver the inner eigensolver to be used.

Member Function Documentation

◆ Solve() [1/2]

int chrono::modal::ChModalSolverDamped::Solve ( const ChAssembly assembly,
ChMatrixDynamic< std::complex< double >> &  eigvects,
ChVectorDynamic< std::complex< double >> &  eigvals,
ChVectorDynamic< double > &  freq,
ChVectorDynamic< double > &  damp_ratios 
) const
virtual

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/2]

int chrono::modal::ChModalSolverDamped::Solve ( const ChSparseMatrix K,
const ChSparseMatrix R,
const ChSparseMatrix M,
const ChSparseMatrix Cq,
ChMatrixDynamic< std::complex< double >> &  eigvects,
ChVectorDynamic< std::complex< double >> &  eigvals,
ChVectorDynamic< double > &  freq,
ChVectorDynamic< double > &  damp_ratios 
) const
virtual

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.


The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono_modal/ChModalSolverDamped.h
  • /builds/uwsbel/chrono/src/chrono_modal/ChModalSolverDamped.cpp