libMesh
adaptive_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_ADAPTIVE_TIME_SOLVER_H
21 #define LIBMESH_ADAPTIVE_TIME_SOLVER_H
22 
23 // Local includes
24 #include "libmesh/system_norm.h"
25 #include "libmesh/first_order_unsteady_solver.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
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::AdaptiveTimeSolver::solve
virtual void solve() override=0
This method solves for the solution at the next timestep.
libMesh::FirstOrderUnsteadySolver
Generic class from which first order UnsteadySolvers should subclass.
Definition: first_order_unsteady_solver.h:74
libMesh::AdaptiveTimeSolver::calculate_norm
virtual Real calculate_norm(System &, NumericVector< Number > &)
A helper function to calculate error norms.
Definition: adaptive_time_solver.C:155
libMesh::AdaptiveTimeSolver::nonlocal_residual
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
Definition: adaptive_time_solver.C:131
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::AdaptiveTimeSolver::AdaptiveTimeSolver
AdaptiveTimeSolver(sys_type &s)
Constructor.
Definition: adaptive_time_solver.C:27
libMesh::AdaptiveTimeSolver::upper_tolerance
Real upper_tolerance
This tolerance is the maximum relative error between an exact time integration and a single time step...
Definition: adaptive_time_solver.h:161
libMesh::DifferentiableSystem
This class provides a specific system class.
Definition: diff_system.h:53
libMesh::AdaptiveTimeSolver::linear_solver
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver() override
An implicit linear solver to use for adjoint and sensitivity problems.
Definition: adaptive_time_solver.C:148
libMesh::AdaptiveTimeSolver::advance_timestep
virtual void advance_timestep() override
This method advances the solution to the next timestep, after a solve() has been performed.
Definition: adaptive_time_solver.C:87
libMesh::AdaptiveTimeSolver::Parent
FirstOrderUnsteadySolver Parent
The parent class.
Definition: adaptive_time_solver.h:55
libMesh::NumericVector< Number >
libMesh::AdaptiveTimeSolver::core_time_solver
std::unique_ptr< UnsteadySolver > core_time_solver
This object is used to take timesteps.
Definition: adaptive_time_solver.h:114
libMesh::AdaptiveTimeSolver::component_norm
SystemNorm component_norm
Error calculations are done in this norm, DISCRETE_L2 by default.
Definition: adaptive_time_solver.h:119
libMesh::AdaptiveTimeSolver::last_deltat
Real last_deltat
We need to store the value of the last deltat used, so that advance_timestep() will increment the sys...
Definition: adaptive_time_solver.h:207
libMesh::AdaptiveTimeSolver::global_tolerance
bool global_tolerance
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accura...
Definition: adaptive_time_solver.h:198
libMesh::AdaptiveTimeSolver::max_growth
Real max_growth
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old de...
Definition: adaptive_time_solver.h:181
libMesh::AdaptiveTimeSolver::error_order
virtual Real error_order() const override
This method is passed on to the core_time_solver.
Definition: adaptive_time_solver.C:102
libMesh::AdaptiveTimeSolver::side_residual
virtual bool side_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
Definition: adaptive_time_solver.C:121
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::AdaptiveTimeSolver
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with ...
Definition: adaptive_time_solver.h:49
libMesh::AdaptiveTimeSolver::max_deltat
Real max_deltat
Do not allow the adaptive time solver to select deltat > max_deltat.
Definition: adaptive_time_solver.h:167
libMesh::AdaptiveTimeSolver::component_scale
std::vector< float > component_scale
If component_norms is non-empty, each variable's contribution to the error of a system will also be s...
Definition: adaptive_time_solver.h:127
libMesh::AdaptiveTimeSolver::element_residual
virtual bool element_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
Definition: adaptive_time_solver.C:111
libMesh::AdaptiveTimeSolver::init
virtual void init() override
The initialization function.
Definition: adaptive_time_solver.C: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::AdaptiveTimeSolver::min_deltat
Real min_deltat
Do not allow the adaptive time solver to select deltat < min_deltat.
Definition: adaptive_time_solver.h:173
libMesh::AdaptiveTimeSolver::reinit
virtual void reinit() override
The reinitialization function.
Definition: adaptive_time_solver.C:77
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::AdaptiveTimeSolver::diff_solver
virtual std::unique_ptr< DiffSolver > & diff_solver() override
An implicit linear or nonlinear solver to use at each timestep.
Definition: adaptive_time_solver.C:141
libMesh::AdaptiveTimeSolver::target_tolerance
Real target_tolerance
This tolerance is the target relative error between an exact time integration and a single time step ...
Definition: adaptive_time_solver.h:144
libMesh::AdaptiveTimeSolver::~AdaptiveTimeSolver
virtual ~AdaptiveTimeSolver()
Destructor.
Definition: adaptive_time_solver.C:47