47 #include <solution/system_of_eqn/SystemOfEqn.h> 48 #include "/usr/include/boost/numeric/ublas/matrix_sparse.hpp" 66 typedef boost::numeric::ublas::mapped_matrix<double> sparse_matrix;
76 void resize_mass_matrix_if_needed(
const size_t &);
84 virtual int solve(
int numModes);
85 virtual int solve(
void);
88 virtual int addA(
const Matrix &,
const ID &,
double fact = 1.0) = 0;
89 virtual int addM(
const Matrix &,
const ID &,
double fact = 1.0) = 0;
91 virtual int setSize(
Graph &theGraph) = 0;
92 virtual void zeroA(
void) = 0;
93 virtual void zeroM(
void);
const int & getNumModes(void) const
Returns the number of computed eigenvalues.
Definition: EigenSOE.cpp:235
virtual const Vector & getEigenvector(int mode) const
Return the autovector that correspond to the mode being passed as parameter.
Definition: EigenSOE.cpp:178
int getNumEqn(void) const
Returns the number of equations.
Definition: EigenSOE.cpp:111
Float vector abstraction.
Definition: Vector.h:93
Matrix getEigenvectors(void) const
Returns a matrix with the computed eigenvectors disposed by columns.
Definition: EigenSOE.cpp:190
virtual void zeroM(void)
Anula la matrix M.
Definition: EigenSOE.cpp:159
Solution procedure for the finite element problem.
Definition: AnalysisAggregation.h:89
double getAngularFrequency(int mode) const
Returns the angular frequency of the i-th mode.
Definition: EigenSOE.cpp:203
Base class for eigenproblem systems of equations.
Definition: EigenSOE.h:63
virtual bool setSolver(EigenSolver *)
Set the solver.
Definition: EigenSOE.cpp:89
Vector getNormalizedEigenvector(int mode) const
Returns the normalized autovector that correspond to the mode being passed as parameter.
Definition: EigenSOE.cpp:185
Matrix getNormalizedEigenvectors(void) const
Returns a matrix whit the normalized eigenvectors disposed by columns (infinity norm).
Definition: EigenSOE.cpp:195
int size
order of A
Definition: EigenSOE.h:68
Matrix getDistributionFactors(void) const
Returns a matrix with the computed distribution factors placed by columns.
Definition: EigenSOE.cpp:274
Vector of integers.
Definition: ID.h:93
Eigenvalue SOE solver.
Definition: EigenSolver.h:59
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies for each mode.
Definition: EigenSOE.cpp:220
virtual ~EigenSOE(void)
Destructor.
Definition: EigenSOE.cpp:134
virtual int solve(void)
No hace nada.
Definition: EigenSOE.cpp:148
Vector getEquivalentStaticLoad(int mode, const double &) const
Return the equivalent static force for the mode passed as parameter.
Definition: EigenSOE.cpp:330
double getEffectiveModalMass(int mode) const
Return the effective modal mass for the i-th mode.
Definition: EigenSOE.cpp:296
double getFrecuencia(int mode) const
Return the frecuency of the i-th mode.
Definition: EigenSOE.cpp:211
System of equations base class.
Definition: SystemOfEqn.h:89
Vector getEigenvalues(void) const
Returns a vector with computed eigenvalues for each mode.
Definition: EigenSOE.cpp:215
Vector getEffectiveModalMasses(void) const
Returns the effective modal masses for each mode.
Definition: EigenSOE.cpp:310
Vector getDistributionFactor(int mode) const
Returns the distribution factors for the i-th mode.
Definition: EigenSOE.cpp:269
double getTotalMass(void) const
Return the model total mass.
Definition: EigenSOE.cpp:320
The Graph class provides the abstraction of a graph.
Definition: Graph.h:93
virtual const double & getEigenvalue(int mode) const
Returns the eigenvalue of the mode passed as parameter.
Definition: EigenSOE.cpp:199
EigenSOE(AnalysisAggregation *, int classTag)
Constructor.
Definition: EigenSOE.cpp:62
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
Matrix of floats.
Definition: Matrix.h:108
double getPeriodo(int mode) const
Returns the period of the i-th mode.
Definition: EigenSOE.cpp:207
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor for the mode.
Definition: EigenSOE.cpp:239
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: EigenSOE.cpp:259
virtual void identityM(void)
Makes M the identity matrix (to find stiffness matrix eigenvalues).
Definition: EigenSOE.cpp:165
Vector getFrecuencias(void) const
Returns a vector with the computed frequencies for each mode.
Definition: EigenSOE.cpp:231
Vector getPeriodos(void) const
Returns a vector with the computed periods for each mode.
Definition: EigenSOE.cpp:226
EigenSolver * getSolver(void)
Return a pointer to the solver used to solve the eigenproblem.
Definition: EigenSOE.cpp:155
sparse_matrix massMatrix
Mass matrix (used in getModalParticipationFactor).
Definition: EigenSOE.h:70