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