adaptive_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_ADAPTIVE_TIME_SOLVER_H
21 #define LIBMESH_ADAPTIVE_TIME_SOLVER_H
22 
23 // Local includes
24 #include "libmesh/system_norm.h"
26 
27 // C++ includes
28 
29 namespace libMesh
30 {
31 
32 // Forward declarations
33 class System;
34 
50 {
51 public:
56 
61  explicit
63 
67  virtual ~AdaptiveTimeSolver ();
68 
69  virtual void init() override;
70 
71  virtual void reinit() override;
72 
73  virtual void solve() override = 0;
74 
75  virtual void advance_timestep() override;
76 
80  virtual Real error_order () const override;
81 
85  virtual bool element_residual (bool get_jacobian,
86  DiffContext &) override;
87 
91  virtual bool side_residual (bool get_jacobian,
92  DiffContext &) override;
93 
97  virtual bool nonlocal_residual (bool get_jacobian,
98  DiffContext &) override;
99 
103  virtual std::unique_ptr<DiffSolver> & diff_solver() override;
104 
109  virtual std::unique_ptr<LinearSolver<Number>> & linear_solver() override;
110 
114  std::unique_ptr<UnsteadySolver> core_time_solver;
115 
120 
127  std::vector<float> component_scale;
128 
145 
162 
168 
174 
182 
199 
200 protected:
201 
208 
213 };
214 
215 
216 } // namespace libMesh
217 
218 
219 #endif // LIBMESH_ADAPTIVE_TIME_SOLVER_H
virtual Real calculate_norm(System &, NumericVector< Number > &)
std::vector< float > component_scale
virtual bool side_residual(bool get_jacobian, DiffContext &) override
virtual Real error_order() const override
virtual std::unique_ptr< DiffSolver > & diff_solver() override
FirstOrderUnsteadySolver Parent
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void reinit() override
virtual void advance_timestep() override
virtual bool element_residual(bool get_jacobian, DiffContext &) override
virtual void solve() override=0
std::unique_ptr< UnsteadySolver > core_time_solver
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver() override
virtual void init() override