66 #ifndef EigenAnalysis_h 67 #define EigenAnalysis_h 69 #include <solution/analysis/analysis/Analysis.h> 87 virtual int analyze(
int numModes);
Float vector abstraction.
Definition: Vector.h:93
EigenAnalysis(AnalysisAggregation *analysis_aggregation)
Constructor.
Definition: EigenAnalysis.cpp:82
double getTotalMass(void) const
Return the masa total del modelo.
Definition: EigenAnalysis.cpp:409
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor corresponding to i mode.
Definition: EigenAnalysis.cpp:341
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies for each mode.
Definition: EigenAnalysis.cpp:301
Solution procedure for the finite element problem.
Definition: AnalysisAggregation.h:89
Base class for eigenproblem systems of equations.
Definition: EigenSOE.h:63
virtual int setIntegrator(EigenIntegrator &theIntegrator)
Sets the integrator to use in the analysis.
Definition: EigenAnalysis.cpp:205
double getAngularFrequency(int mode) const
Return the angular frequency for the i-th mode.
Definition: EigenAnalysis.cpp:279
Matrix getNormalizedEigenvectors(void) const
Returns a matrix con los eigenvectors normalizados (infinity norm) as columns of the matrix...
Definition: EigenAnalysis.cpp:258
void clearAll(void)
Clears all object members (constraint handler, analysis model,...).
Definition: EigenAnalysis.cpp:90
virtual int analyze(int numModes)
Performs the analysis.
Definition: EigenAnalysis.cpp:100
Base class for the object that perform the analysis.
Definition: Analysis.h:116
Vector getEquivalentStaticLoad(int mode, const double &) const
Returns the equivalent static load for the mode being passed as parameter.
Definition: EigenAnalysis.cpp:420
Eigenproblem analysis.
Definition: EigenAnalysis.h:78
virtual int setEigenSOE(EigenSOE &theSOE)
Sets the sistema de eigenvalues to use in the analysis.
Definition: EigenAnalysis.cpp:214
double getPeriodo(int mode) const
Returns the period for the i-th mode.
Definition: EigenAnalysis.cpp:283
virtual const Vector & getEigenvector(int mode) const
Returns the autovector that corresponds to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:222
Vector getEffectiveModalMasses(void) const
Returns the masas modales efectivas.
Definition: EigenAnalysis.cpp:399
Vector getNormalizedEigenvector(int mode) const
Returns the normalized autovector that correspond to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:235
virtual const double & getEigenvalue(int mode) const
Returns the eigenvalue that corresponds to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:269
Solution algorithm for eigenproblem.
Definition: EigenAlgorithm.h:83
virtual int domainChanged(void)
Hace los cambios necesarios como consecuencia de un cambio en el domain.
Definition: EigenAnalysis.cpp:141
Base class for eigenproblem integrators.
Definition: EigenIntegrator.h:85
Matrix getEigenvectors(void) const
Returns a matrix with the computed eigenvectors as columns of the matrix.
Definition: EigenAnalysis.cpp:247
Vector getDistributionFactor(int mode) const
Returns the distribution factor correspondint to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:367
int getNumModes(void) const
Returns the number of eigenvalues que se han calculado.
Definition: EigenAnalysis.cpp:331
Solution procedure for the finite element problem.
Definition: ProcSolu.h:56
Vector getFrecuencias(void) const
Returns a vector con las frecuencias calculadas.
Definition: EigenAnalysis.cpp:321
double getEffectiveModalMass(int mode) const
Return the masa modal efectiva correspondiente al modo i.
Definition: EigenAnalysis.cpp:389
Vector getEigenvalues(void) const
Returns a vector with the computed eigenvalues for each mode.
Definition: EigenAnalysis.cpp:291
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: EigenAnalysis.cpp:354
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:34
Matrix of floats.
Definition: Matrix.h:108
Matrix getDistributionFactors(void) const
Returns a matrix with the distribution factors by columns.
Definition: EigenAnalysis.cpp:378
virtual int setAlgorithm(EigenAlgorithm &theAlgo)
Sets the algorithm to use in the analysis.
Definition: EigenAnalysis.cpp:196
Analysis * getCopy(void) const
Virtual constructor.
Definition: EigenAnalysis.cpp:86
double getFrecuencia(int mode) const
Return the frequency for the i-th mode.
Definition: EigenAnalysis.cpp:287
Vector getPeriodos(void) const
Returns a vector with the computed vectors for each mode.
Definition: EigenAnalysis.cpp:311