eigen_time_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_EIGEN_TIME_SOLVER_H
21 #define LIBMESH_EIGEN_TIME_SOLVER_H
22 
23 #include "libmesh/libmesh_config.h"
24 #ifdef LIBMESH_HAVE_SLEPC
25 
26 // Local includes
27 #include "libmesh/time_solver.h"
28 
29 // C++ includes
30 
31 namespace libMesh
32 {
33 
34 // Forward declarations
35 template <typename T> class EigenSolver;
36 
37 
67 {
68 public:
73 
77  typedef TimeSolver Parent;
78 
83  explicit
85 
89  virtual ~EigenTimeSolver ();
90 
95  virtual void init () override;
96 
101  virtual void reinit () override;
102 
107  virtual void solve () override;
108 
112  virtual void advance_timestep () override {}
113 
118  Real error_order() const { return 0.; }
119 
124  virtual bool element_residual (bool get_jacobian,
125  DiffContext &) override;
126 
130  virtual bool side_residual (bool get_jacobian,
131  DiffContext &) override;
132 
136  virtual bool nonlocal_residual (bool get_jacobian,
137  DiffContext &) override;
138 
144  virtual Real du (const SystemNorm &) const override { return 0.; }
145 
149  virtual bool is_steady() const override { return true; }
150 
155  std::unique_ptr<EigenSolver<Number>> eigen_solver;
156 
162 
166  unsigned int maxits;
167 
172 
183 
189 
194  unsigned int n_iterations_reqd;
195 
196 private:
197 
203 
208 
213  };
214 
219 };
220 
221 } // namespace libMesh
222 
223 
224 #endif // LIBMESH_HAVE_SLEPC
225 #endif // LIBMESH_EIGEN_TIME_SOLVER_H
virtual void init() override
virtual void solve() override
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
virtual Real du(const SystemNorm &) const override
virtual bool is_steady() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DifferentiableSystem sys_type
virtual bool side_residual(bool get_jacobian, DiffContext &) override
virtual void reinit() override
virtual bool element_residual(bool get_jacobian, DiffContext &) override
unsigned int n_eigenpairs_to_compute
std::unique_ptr< EigenSolver< Number > > eigen_solver
virtual void advance_timestep() override