libMesh
time_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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_TIME_SOLVER_H
21 #define LIBMESH_TIME_SOLVER_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/linear_solver.h"
26 #include "libmesh/numeric_vector.h"
27 #include "libmesh/reference_counted_object.h"
28 #include "libmesh/solution_history.h"
29 
30 // C++ includes
31 #include <memory>
32 
33 namespace libMesh
34 {
35 
36 // Forward Declarations
37 class DiffContext;
38 class DiffSolver;
39 class DifferentiablePhysics;
40 class DifferentiableSystem;
41 class ParameterVector;
42 class SystemNorm;
43 
58 class TimeSolver : public ReferenceCountedObject<TimeSolver>
59 {
60 public:
65 
70  explicit
71  TimeSolver (sys_type & s);
72 
76  virtual ~TimeSolver ();
77 
82  virtual void init ();
83 
89  virtual void init_data ();
90 
95  virtual void reinit ();
96 
103  virtual void solve ();
104 
111  virtual void advance_timestep ();
112 
118  virtual void adjoint_advance_timestep ();
119 
124  virtual void retrieve_timestep();
125 
135  virtual bool element_residual (bool request_jacobian,
136  DiffContext &) = 0;
137 
147  virtual bool side_residual (bool request_jacobian,
148  DiffContext &) = 0;
149 
159  virtual bool nonlocal_residual (bool request_jacobian,
160  DiffContext &) = 0;
161 
166  virtual void before_timestep () {}
167 
171  const sys_type & system () const { return _system; }
172 
176  sys_type & system () { return _system; }
177 
181  virtual std::unique_ptr<DiffSolver> & diff_solver() { return _diff_solver; }
182 
186  virtual std::unique_ptr<LinearSolver<Number>> & linear_solver() { return _linear_solver; }
187 
191  bool quiet;
192 
202  virtual Real du(const SystemNorm & norm) const = 0;
203 
207  virtual bool is_steady() const = 0;
208 
221 
226  void set_solution_history(const SolutionHistory & _solution_history);
227 
232  bool is_adjoint() const
233  { return _is_adjoint; }
234 
239  void set_is_adjoint(bool _is_adjoint_value)
240  { _is_adjoint = _is_adjoint_value; }
241 
242 protected:
243 
247  std::unique_ptr<DiffSolver> _diff_solver;
248 
252  std::unique_ptr<LinearSolver<Number>> _linear_solver;
253 
258 
264  std::unique_ptr<SolutionHistory> solution_history;
265 
270  typedef bool (DifferentiablePhysics::*ResFuncType) (bool, DiffContext &);
271 
272  typedef void (DiffContext::*ReinitFuncType) (Real);
273 
274 private:
275 
281 
282 };
283 
284 
285 } // namespace libMesh
286 
287 
288 #endif // LIBMESH_TIME_SOLVER_H
libMesh::TimeSolver::sys_type
DifferentiableSystem sys_type
The type of system.
Definition: time_solver.h:64
libMesh::TimeSolver::reinit
virtual void reinit()
The reinitialization function.
Definition: time_solver.C:47
libMesh::TimeSolver::ResFuncType
bool(DifferentiablePhysics::* ResFuncType)(bool, DiffContext &)
Definitions of argument types for use in refactoring subclasses.
Definition: time_solver.h:270
libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure
unsigned int reduce_deltat_on_diffsolver_failure
This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat ...
Definition: time_solver.h:220
libMesh::TimeSolver::~TimeSolver
virtual ~TimeSolver()
Destructor.
Definition: time_solver.C:41
libMesh::TimeSolver::init
virtual void init()
The initialization function.
Definition: time_solver.C:63
libMesh::DifferentiablePhysics
This class provides a specific system class.
Definition: diff_physics.h:76
libMesh::TimeSolver::quiet
bool quiet
Print extra debugging information if quiet == false.
Definition: time_solver.h:191
libMesh::TimeSolver::_diff_solver
std::unique_ptr< DiffSolver > _diff_solver
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:247
libMesh::ReferenceCountedObject
This class implements reference counting.
Definition: reference_counted_object.h:65
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::TimeSolver::system
sys_type & system()
Definition: time_solver.h:176
libMesh::TimeSolver::ReinitFuncType
void(DiffContext::* ReinitFuncType)(Real)
Definition: time_solver.h:272
libMesh::TimeSolver::advance_timestep
virtual void advance_timestep()
This method advances the solution to the next timestep, after a solve() has been performed.
Definition: time_solver.C:101
libMesh::DifferentiableSystem
This class provides a specific system class.
Definition: diff_system.h:53
libMesh::TimeSolver::element_residual
virtual bool element_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics element_time_derivative(), element_constraint(),...
libMesh::TimeSolver::adjoint_advance_timestep
virtual void adjoint_advance_timestep()
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been...
Definition: time_solver.C:105
libMesh::TimeSolver::set_solution_history
void set_solution_history(const SolutionHistory &_solution_history)
A setter function users will employ if they need to do something other than save no solution history.
Definition: time_solver.C:96
libMesh::TimeSolver::_is_adjoint
bool _is_adjoint
This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.
Definition: time_solver.h:280
libMesh::TimeSolver::TimeSolver
TimeSolver(sys_type &s)
Constructor.
Definition: time_solver.C:28
libMesh::TimeSolver::linear_solver
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
An implicit linear solver to use for adjoint and sensitivity problems.
Definition: time_solver.h:186
libMesh::TimeSolver::nonlocal_residual
virtual bool nonlocal_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics nonlocal_time_derivative(), nonlocal_constraint(),...
libMesh::TimeSolver::diff_solver
virtual std::unique_ptr< DiffSolver > & diff_solver()
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:181
libMesh::TimeSolver::set_is_adjoint
void set_is_adjoint(bool _is_adjoint_value)
Accessor for setting whether we need to do a primal or adjoint solve.
Definition: time_solver.h:239
libMesh::TimeSolver::_system
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:257
libMesh::TimeSolver::before_timestep
virtual void before_timestep()
This method is for subclasses or users to override to do arbitrary processing between timesteps.
Definition: time_solver.h:166
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::TimeSolver::init_data
virtual void init_data()
The data initialization function.
Definition: time_solver.C:76
libMesh::TimeSolver::is_adjoint
bool is_adjoint() const
Accessor for querying whether we need to do a primal or adjoint solve.
Definition: time_solver.h:232
libMesh::TimeSolver::du
virtual Real du(const SystemNorm &norm) const =0
Computes the size of ||u^{n+1} - u^{n}|| in some norm.
libMesh::TimeSolver::retrieve_timestep
virtual void retrieve_timestep()
This method retrieves all the stored solutions at the current system.time.
Definition: time_solver.C:109
std::norm
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::TimeSolver
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:58
libMesh::SystemNorm
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:51
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::TimeSolver::system
const sys_type & system() const
Definition: time_solver.h:171
libMesh::TimeSolver::side_residual
virtual bool side_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics side_time_derivative(), side_constraint(),...
libMesh::TimeSolver::solution_history
std::unique_ptr< SolutionHistory > solution_history
A std::unique_ptr to a SolutionHistory object.
Definition: time_solver.h:264
libMesh::SolutionHistory
A SolutionHistory class that enables the storage and retrieval of timesteps and (in the future) adapt...
Definition: solution_history.h:35
libMesh::TimeSolver::is_steady
virtual bool is_steady() const =0
Is this effectively a steady-state solver?
libMesh::TimeSolver::_linear_solver
std::unique_ptr< LinearSolver< Number > > _linear_solver
An implicit linear solver to use for adjoint problems.
Definition: time_solver.h:252
libMesh::TimeSolver::solve
virtual void solve()
This method solves for the solution at the next timestep (or solves for a steady-state solution).
Definition: time_solver.C:88