xc
EigenSOE.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 // File: ~/system_of_eqn/eigenSOE/EigenSOE.h
29 //
30 // Written: Jun Peng
31 // Created: Sat Feb. 6, 1999
32 // Revision: A
33 //
34 // Description: This file contains the class definition of EigenSOE.
35 // EigenSOE is a subclass of SystemOfEqn.
36 // It has pure virtual functions which must be implemented in it's derived
37 // subclasses. To solve the general eigen value equations means that
38 // by the given K and M, find the corresponding eigen value and eigen
39 // vectors.
40 //
41 // This class is inherited from the base class of SystemOfEqn
42 // which was created by fmk (Frank).
43 
44 
45 #ifndef EigenSOE_h
46 #define EigenSOE_h
47 
48 #include <solution/system_of_eqn/SystemOfEqn.h>
49 #include "/usr/include/boost/numeric/ublas/matrix_sparse.hpp"
50 
51 namespace XC {
52 class EigenSolver;
53 class Matrix;
54 class Vector;
55 class ID;
56 
60 //
62 //
64 class EigenSOE: public SystemOfEqn
65  {
66  public:
67  typedef boost::numeric::ublas::mapped_matrix<double> sparse_matrix;
68  protected:
69  int size;
70  bool factored;
71  sparse_matrix massMatrix;
72  EigenSolver *theSolver;
73 
74  void free_mem(void);
75  void copy(const EigenSolver *);
76  virtual bool setSolver(EigenSolver *);
77  void resize_mass_matrix_if_needed(const size_t &);
78 
79  EigenSOE(SolutionStrategy *,int classTag);
80  public:
81  virtual ~EigenSOE(void);
82 
83  int getNumEqn(void) const;
84 
85  virtual int solve(int numModes);
86  virtual int solve(void);
87  virtual double getDeterminant(void);
88  virtual double getRCond(const char &norm= '1');
89 
90  // pure virtual functions
91  virtual int addA(const Matrix &, const ID &, double fact = 1.0) = 0;
92  virtual int addM(const Matrix &, const ID &, double fact = 1.0) = 0;
93 
94  virtual int setSize(Graph &theGraph) = 0;
95  virtual void zeroA(void) = 0;
96  virtual void zeroM(void);
97  virtual void identityM(void);
98  virtual boost::python::list getMPy(void) const;
99 
100  const int &getNumModes(void) const;
101  bool standardProblem(void) const;
102 
103  //Eigenvectors
104  virtual const Vector &getEigenvector(int mode) const;
105  Vector getNormalizedEigenvector(int mode) const;
106  Matrix getEigenvectors(void) const;
107  Matrix getNormalizedEigenvectors(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(const int &numModes) const;
115  Vector getEigenvalues(void) const;
116  Vector getAngularFrequencies(void) const;
117  Vector getPeriods(void) const;
118  Vector getFrequencies(void) const;
119 
120  //Modal participation factors.
121  virtual double getModalParticipationFactor(int mode) const;
122  Vector getModalParticipationFactors(const int &numModes) 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  EigenSolver *getSolver(void);
138  EigenSolver &newSolver(const std::string &);
139  };
140 } // end of XC namespace
141 
142 #endif
143 
144 
145 
const int & getNumModes(void) const
Returns the number of computed eigenvalues.
Definition: EigenSOE.cpp:295
virtual const Vector & getEigenvector(int mode) const
Return the autovector that correspond to the mode being passed as parameter.
Definition: EigenSOE.cpp:229
int getNumEqn(void) const
Returns the number of equations.
Definition: EigenSOE.cpp:116
EigenSOE(SolutionStrategy *, int classTag)
Constructor.
Definition: EigenSOE.cpp:65
Float vector abstraction.
Definition: Vector.h:94
virtual boost::python::list getMPy(void) const
Return the rows of the mass matrix in a Python list.
Definition: EigenSOE.cpp:194
Matrix getEigenvectors(void) const
Returns a matrix with the computed eigenvectors disposed by columns.
Definition: EigenSOE.cpp:241
Vector getFrequencies(void) const
Returns a vector with the computed frequencies for each mode.
Definition: EigenSOE.cpp:291
virtual void zeroM(void)
Anula la matrix M.
Definition: EigenSOE.cpp:188
double getAngularFrequency(int mode) const
Returns the angular frequency of the i-th mode.
Definition: EigenSOE.cpp:254
Base class for eigenproblem systems of equations.
Definition: EigenSOE.h:64
virtual bool setSolver(EigenSolver *)
Set the solver.
Definition: EigenSOE.cpp:93
Vector getNormalizedEigenvector(int mode) const
Returns the normalized autovector that correspond to the mode being passed as parameter.
Definition: EigenSOE.cpp:236
Matrix getNormalizedEigenvectors(void) const
Returns a matrix whit the normalized eigenvectors disposed by columns (infinity norm).
Definition: EigenSOE.cpp:246
int size
order of A
Definition: EigenSOE.h:69
Matrix getDistributionFactors(void) const
Returns a matrix with the computed distribution factors placed by columns.
Definition: EigenSOE.cpp:352
virtual double getDeterminant(void)
Returns the determinant of the system matrix.
Definition: EigenSOE.cpp:168
Vector of integers.
Definition: ID.h:95
Eigenvalue SOE solver.
Definition: EigenSolver.h:60
Vector getAngularFrequencies(void) const
Returns a vector with the computed angular frequencies for each mode.
Definition: EigenSOE.cpp:280
virtual ~EigenSOE(void)
Destructor.
Definition: EigenSOE.cpp:146
virtual int solve(void)
Do nothing.
Definition: EigenSOE.cpp:160
virtual double getRCond(const char &norm='1')
Returns the reciprocal of the condition number.
Definition: EigenSOE.cpp:178
Vector getEquivalentStaticLoad(int mode, const double &) const
Return the equivalent static force for the mode passed as parameter.
Definition: EigenSOE.cpp:408
double getEffectiveModalMass(int mode) const
Return the effective modal mass for the i-th mode.
Definition: EigenSOE.cpp:374
System of equations base class.
Definition: SystemOfEqn.h:90
Vector getEigenvalues(void) const
Returns a vector with computed eigenvalues for each mode.
Definition: EigenSOE.cpp:275
Vector getEffectiveModalMasses(void) const
Returns the effective modal masses for each mode.
Definition: EigenSOE.cpp:388
Vector getDistributionFactor(int mode) const
Returns the distribution factors for the i-th mode.
Definition: EigenSOE.cpp:347
double getTotalMass(void) const
Return the model total mass.
Definition: EigenSOE.cpp:398
The Graph class provides the abstraction of a graph.
Definition: Graph.h:94
virtual const double & getEigenvalue(int mode) const
Returns the eigenvalue of the mode passed as parameter.
Definition: EigenSOE.cpp:250
Vector getPeriods(void) const
Returns a vector with the computed periods for each mode.
Definition: EigenSOE.cpp:286
Solution strategy for the finite element problem.
Definition: SolutionStrategy.h:94
double getFrequency(int mode) const
Return the frequency of the i-th mode.
Definition: EigenSOE.cpp:262
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
virtual double getModalParticipationFactor(int mode) const
Returns the modal participation factor for the mode.
Definition: EigenSOE.cpp:308
Vector getModalParticipationFactors(void) const
Returns the modal participation factors.
Definition: EigenSOE.cpp:340
EigenSolver & newSolver(const std::string &)
Set the solver to use.
Definition: EigenSOE.cpp:121
virtual void identityM(void)
Makes M the identity matrix (to find stiffness matrix eigenvalues).
Definition: EigenSOE.cpp:216
EigenSolver * getSolver(void)
Return a pointer to the solver used to solve the eigenproblem.
Definition: EigenSOE.cpp:183
double getPeriod(int mode) const
Returns the period of the i-th mode.
Definition: EigenSOE.cpp:258
sparse_matrix massMatrix
Mass matrix (used in getModalParticipationFactor).
Definition: EigenSOE.h:71