67 #ifndef EigenAnalysis_h 68 #define EigenAnalysis_h 70 #include <solution/analysis/analysis/Analysis.h> 88 virtual int analyze(
int numModes);
96 int getDomainStamp(
void)
const 97 {
return domainStamp; }
98 void setDomainStamp(
const int &i)
Float vector abstraction.
Definition: Vector.h:94
double getPeriod(int mode) const
Returns the period for the i-th mode.
Definition: EigenAnalysis.cpp:306
double getTotalMass(void) const
Return the masa total del modelo.
Definition: EigenAnalysis.cpp:439
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor corresponding to i mode.
Definition: EigenAnalysis.cpp:371
boost::python::list getNormalizedEigenvectorsPy(void) const
Returns a Python list with the computed eigenvectors as lists.
Definition: EigenAnalysis.cpp:259
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies for each mode.
Definition: EigenAnalysis.cpp:331
Base class for eigenproblem systems of equations.
Definition: EigenSOE.h:64
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:302
Matrix getNormalizedEigenvectors(void) const
Returns a matrix with the normalized eigenvectors (infinity norm) as columns of the matrix...
Definition: EigenAnalysis.cpp:281
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:117
Vector getEquivalentStaticLoad(int mode, const double &) const
Returns the equivalent static load for the mode being passed as parameter.
Definition: EigenAnalysis.cpp:450
Eigenproblem analysis.
Definition: EigenAnalysis.h:79
Vector getFrequencies(void) const
Returns a vector with the computed frequencies.
Definition: EigenAnalysis.cpp:351
virtual int setEigenSOE(EigenSOE &theSOE)
Sets the system of eigenvalues to use in the analysis.
Definition: EigenAnalysis.cpp:216
virtual const Vector & getEigenvector(int mode) const
Returns the autovector that corresponds to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:224
EigenAnalysis(SolutionStrategy *analysis_aggregation)
Constructor.
Definition: EigenAnalysis.cpp:82
Vector getEffectiveModalMasses(void) const
Returns the masas modales efectivas.
Definition: EigenAnalysis.cpp:429
Solution procedure for the finite element problem.
Definition: SolutionProcedure.h:57
Vector getNormalizedEigenvector(int mode) const
Returns the normalized autovector that correspond to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:237
virtual const double & getEigenvalue(int mode) const
Returns the eigenvalue that corresponds to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:292
Solution algorithm for eigenproblem.
Definition: EigenAlgorithm.h:84
virtual int domainChanged(void)
Make the changes derived of a change in the domain.
Definition: EigenAnalysis.cpp:141
double getFrequency(int mode) const
Return the frequency for the i-th mode.
Definition: EigenAnalysis.cpp:310
Base class for eigenproblem integrators.
Definition: EigenIntegrator.h:86
Matrix getEigenvectors(void) const
Returns a matrix with the computed eigenvectors as columns of the matrix.
Definition: EigenAnalysis.cpp:249
Vector getDistributionFactor(int mode) const
Returns the distribution factor correspondint to the mode being passed as parameter.
Definition: EigenAnalysis.cpp:397
int getNumModes(void) const
Returns the number of computed eigenvalues.
Definition: EigenAnalysis.cpp:361
Solution strategy for the finite element problem.
Definition: SolutionStrategy.h:94
double getEffectiveModalMass(int mode) const
Return the masa modal efectiva correspondiente al modo i.
Definition: EigenAnalysis.cpp:419
Vector getEigenvalues(void) const
Returns a vector with the computed eigenvalues for each mode.
Definition: EigenAnalysis.cpp:314
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: EigenAnalysis.cpp:384
Vector getPeriods(void) const
Returns a vector with the computed vectors for each mode.
Definition: EigenAnalysis.cpp:341
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
Matrix getDistributionFactors(void) const
Returns a matrix with the distribution factors by columns.
Definition: EigenAnalysis.cpp:408
virtual int setAlgorithm(EigenAlgorithm &theAlgo)
Sets the algorithm to use in the analysis.
Definition: EigenAnalysis.cpp:196
boost::python::list getEigenvaluesPy(void) const
Returns a Python list with the computed eigenvalues for each mode.
Definition: EigenAnalysis.cpp:324
boost::python::list getNormalizedEigenvectorPy(int mode) const
Return a Python list with the component of the eigenvector corrsponding to the given mode...
Definition: EigenAnalysis.cpp:273
Analysis * getCopy(void) const
Virtual constructor.
Definition: EigenAnalysis.cpp:86