xc
EigenAnalysis.h
1 //----------------------------------------------------------------------------
2 // XC program; finite element analysis code
3 // for structural analysis and design.
4 //
5 // Copyright (C) Luis Claudio Pérez Tato
6 //
7 // This program derives from OpenSees <http://opensees.berkeley.edu>
8 // developed by the «Pacific earthquake engineering research center».
9 //
10 // Except for the restrictions that may arise from the copyright
11 // of the original program (see copyright_opensees.txt)
12 // XC is free software: you can redistribute it and/or modify
13 // it under the terms of the GNU General Public License as published by
14 // the Free Software Foundation, either version 3 of the License, or
15 // (at your option) any later version.
16 //
17 // This software is distributed in the hope that it will be useful, but
18 // WITHOUT ANY WARRANTY; without even the implied warranty of
19 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 // GNU General Public License for more details.
21 //
22 //
23 // You should have received a copy of the GNU General Public License
24 // along with this program.
25 // If not, see <http://www.gnu.org/licenses/>.
26 //----------------------------------------------------------------------------
27 /* ****************************************************************** **
28 ** OpenSees - Open System for Earthquake Engineering Simulation **
29 ** Pacific Earthquake Engineering Research Center **
30 ** **
31 ** **
32 ** (C) Copyright 1999, The Regents of the University of California **
33 ** All Rights Reserved. **
34 ** **
35 ** Commercial use of this program without express permission of the **
36 ** University of California, Berkeley, is strictly prohibited. See **
37 ** file 'COPYRIGHT' in main directory for information on usage and **
38 ** redistribution, and for a DISCLAIMER OF ALL WARRANTIES. **
39 ** **
40 ** Developed by: **
41 ** Frank McKenna (fmckenna@ce.berkeley.edu) **
42 ** Gregory L. Fenves (fenves@ce.berkeley.edu) **
43 ** Filip C. Filippou (filippou@ce.berkeley.edu) **
44 ** **
45 ** ****************************************************************** */
46 
47 // $Revision: 1.1.1.1 $
48 // $Date: 2000/09/15 08:23:16 $
49 // $Source: /usr/local/cvs/OpenSees/SRC/analysis/analysis/EigenAnalysis.h,v $
50 
51 
52 // File: ~/analysis/analysis/eigenAnalysis/EigenAnalysis.h
53 //
54 // Written: Jun Peng
55 // Created: Wed Jan 27, 1999
56 // Revision: A
57 //
58 // Description: This file contains the class definition of EigenAnalysis.
59 // EigenAnalysis is a subclass of Analysis, it is used to perform the
60 // eigen value analysis on the FE_Model.
61 //
62 // This class is inheritanted from the base class of Analysis
63 // which was created by fmk (Frank).
64 
65 
66 #ifndef EigenAnalysis_h
67 #define EigenAnalysis_h
68 
69 #include <solution/analysis/analysis/Analysis.h>
70 
71 namespace XC {
72  class Vector;
73  class Matrix;
74 
76 //
78 class EigenAnalysis : public Analysis
79  {
80  protected:
81  int domainStamp;
82 
83  friend class ProcSolu;
84  EigenAnalysis(AnalysisAggregation *analysis_aggregation);
85  Analysis *getCopy(void) const;
86  public:
87  virtual int analyze(int numModes);
88  void clearAll(void);
89  virtual int domainChanged(void);
90 
91  virtual int setAlgorithm(EigenAlgorithm &theAlgo);
92  virtual int setIntegrator(EigenIntegrator &theIntegrator);
93  virtual int setEigenSOE(EigenSOE &theSOE);
94 
95  //Eigenvectors.
96  virtual const Vector &getEigenvector(int mode) const;
97  Vector getNormalizedEigenvector(int mode) const;
98  Matrix getEigenvectors(void) const;
100 
101  //Eigenvalues.
102  virtual const double &getEigenvalue(int mode) const;
103  double getAngularFrequency(int mode) const;
104  double getPeriodo(int mode) const;
105  double getFrecuencia(int mode) const;
106  Vector getEigenvalues(void) const;
107  Vector getAngularFrequencies(void) const;
108  Vector getPeriodos(void) const;
109  Vector getFrecuencias(void) const;
110  int getNumModes(void) const;
111 
112  //Modal participation factors.
113  virtual double getModalParticipationFactor(int mode) const;
115 
116  //Distribution factors.
117  Vector getDistributionFactor(int mode) const;
118  Matrix getDistributionFactors(void) const;
119 
120  //Effective modal masses.
121  double getEffectiveModalMass(int mode) const;
122  Vector getEffectiveModalMasses(void) const;
123  double getTotalMass(void) const;
124 
125  //Equivalent static load.
126  Vector getEquivalentStaticLoad(int mode,const double &) const;
127  };
128 
129 } // end of XC namespace
130 
131 #endif
132 
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