xc
SpectraSOE.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 //SpectraSOE.h
29 
30 #ifndef SpectraSOE_h
31 #define SpectraSOE_h
32 
33 #include <solution/system_of_eqn/eigenSOE/EigenSOE.h>
34 #include <Eigen/Core>
35 #include <Eigen/SparseCore>
36 
37 
38 namespace XC {
39 class SpectraSolver;
40 
42 //
44 class SpectraSOE: public EigenSOE
45  {
46  private:
47  typedef Eigen::Triplet<double> T;
48  std::deque<T> tripletListA;
49  std::deque<T> tripletListM;
50 
51  Eigen::SparseMatrix<double> A;
52  Eigen::SparseMatrix<double> M;
53 
54  int addToMatrix(std::deque<T> &,const Matrix &, const ID &,const double &);
55  void store_mass_matrix(void);
56  protected:
57  bool setSolver(EigenSolver *);
58  void assembleMatrices(void);
59 
60  friend class SolutionStrategy;
61  friend class FEM_ObjectBroker;
63  SystemOfEqn *getCopy(void) const;
64  public:
65  virtual int setSize(Graph &theGraph);
66 
67  virtual int addA(const Matrix &, const ID &, double fact = 1.0);
68  virtual int addM(const Matrix &, const ID &, double fact = 1.0);
69 
70  virtual void zeroA(void);
71  virtual void zeroM(void);
72  virtual void identityM(void);
73 
74  inline const Eigen::SparseMatrix<double> &getA(void) const
75  { return A; }
76  virtual boost::python::list getAPy(void) const;
77  inline const Eigen::SparseMatrix<double> &getM(void) const
78  { return M; }
79  //virtual boost::python::list getMPy(void) const;
80 
81  int sendSelf(Communicator &);
82  int recvSelf(const Communicator &);
83 
84  friend class SpectraSolver;
85  };
86 inline SystemOfEqn *SpectraSOE::getCopy(void) const
87  { return new SpectraSOE(*this); }
88 } // end of XC namespace
89 
90 #endif
Base class for eigenproblem systems of equations.
Definition: EigenSOE.h:64
Communication parameters between processes.
Definition: Communicator.h:66
virtual void identityM(void)
Makes M the identity matrix (to find stiffness matrix eigenvalues).
Definition: SpectraSOE.cc:241
virtual int setSize(Graph &theGraph)
Sets the system size.
Definition: SpectraSOE.cc:56
Arpack++ based band matrix eigenvalue SOE solver.
Definition: SpectraSOE.h:44
bool setSolver(EigenSolver *)
Sets the solver to use.
Definition: SpectraSOE.cc:40
virtual void zeroM(void)
Zeroes the matrix M.
Definition: SpectraSOE.cc:233
FEM_ObjectBroker is is an object broker class for the finite element method.
Definition: FEM_ObjectBroker.h:151
Vector of integers.
Definition: ID.h:95
SpectraSOE(SolutionStrategy *)
Constructor.
Definition: SpectraSOE.cc:36
virtual int addM(const Matrix &, const ID &, double fact=1.0)
Assemblies into M the matrix being passed as parameter multimplied by the fact parameter.
Definition: SpectraSOE.cc:150
Eigenvalue SOE solver.
Definition: EigenSolver.h:60
virtual void zeroA(void)
Zeroes the matrix A.
Definition: SpectraSOE.cc:142
System of equations base class.
Definition: SystemOfEqn.h:90
virtual int addA(const Matrix &, const ID &, double fact=1.0)
Assemblies into A the matrix being passed as parameter multimplied by the fact parameter.
Definition: SpectraSOE.cc:138
The Graph class provides the abstraction of a graph.
Definition: Graph.h:94
Solution strategy for the finite element problem.
Definition: SolutionStrategy.h:94
virtual boost::python::list getAPy(void) const
Return the rows of the mass matrix in a Python list.
Definition: SpectraSOE.cc:269
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
int sendSelf(Communicator &)
Send the object.
Definition: SpectraSOE.cc:290
Arpack++ based band matrix eigenproblem solver.
Definition: SpectraSolver.h:41
Matrix of floats.
Definition: Matrix.h:111
int recvSelf(const Communicator &)
Receive the object.
Definition: SpectraSOE.cc:293
void assembleMatrices(void)
Assemble the matrices.
Definition: SpectraSOE.cc:224
SystemOfEqn * getCopy(void) const
Virtual constructor.
Definition: SpectraSOE.h:86