libMesh
adaptive_time_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 std::pair<unsigned int, Real> adjoint_solve (const QoISet & qoi_indices) override = 0;
76 
77  virtual void advance_timestep() override;
78 
79  virtual void adjoint_advance_timestep() override;
80 
81  virtual void retrieve_timestep() override;
82 
83  virtual void integrate_qoi_timestep() override = 0;
84 
85  virtual void integrate_adjoint_sensitivity(const QoISet & qois, const ParameterVector & parameter_vector, SensitivityData & sensitivities) override = 0;
86 
87 #ifdef LIBMESH_ENABLE_AMR
88  virtual void integrate_adjoint_refinement_error_estimate(AdjointRefinementEstimator & adjoint_refinement_error_estimator, ErrorVector & QoI_elementwise_error) override = 0;
89 #endif // LIBMESH_ENABLE_AMR
90 
92 
96  virtual Real error_order () const override;
97 
101  virtual bool element_residual (bool get_jacobian,
102  DiffContext &) override;
103 
107  virtual bool side_residual (bool get_jacobian,
108  DiffContext &) override;
109 
113  virtual bool nonlocal_residual (bool get_jacobian,
114  DiffContext &) override;
115 
119  virtual std::unique_ptr<DiffSolver> & diff_solver() override;
120 
125  virtual std::unique_ptr<LinearSolver<Number>> & linear_solver() override;
126 
130  std::unique_ptr<UnsteadySolver> core_time_solver;
131 
136 
143  std::vector<float> component_scale;
144 
161 
178 
184 
190 
198 
206 
223 
224 protected:
225 
230 };
231 
232 
233 } // namespace libMesh
234 
235 
236 #endif // LIBMESH_ADAPTIVE_TIME_SOLVER_H
virtual Real calculate_norm(System &, NumericVector< Number > &)
A helper function to calculate error norms.
Real completed_timestep_size
The adaptive time solver&#39;s have two notions of deltat.
std::vector< float > component_scale
If component_norms is non-empty, each variable&#39;s contribution to the error of a system will also be s...
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
virtual bool side_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
SystemNorm component_norm
Error calculations are done in this norm, DISCRETE_L2 by default.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
virtual void integrate_adjoint_sensitivity(const QoISet &qois, const ParameterVector &parameter_vector, SensitivityData &sensitivities) override=0
A method to integrate the adjoint sensitivity w.r.t a given parameter vector.
This class implements a "brute force" goal-oriented error estimator which computes an estimate of err...
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
virtual void adjoint_advance_timestep() override
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been...
virtual Real last_completed_timestep_size() override
Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all subste...
Real target_tolerance
This tolerance is the target relative error between an exact time integration and a single time step ...
virtual Real error_order() const override
This method is passed on to the core_time_solver.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Definition: system_norm.h:45
virtual std::unique_ptr< DiffSolver > & diff_solver() override
An implicit linear or nonlinear solver to use at each timestep.
The libMesh namespace provides an interface to certain functionality in the library.
virtual void integrate_adjoint_refinement_error_estimate(AdjointRefinementEstimator &adjoint_refinement_error_estimator, ErrorVector &QoI_elementwise_error) override=0
A method to compute the adjoint refinement error estimate at the current timestep.
FirstOrderUnsteadySolver Parent
The parent class.
AdaptiveTimeSolver(sys_type &s)
Constructor.
This class provides a specific system class.
Definition: diff_system.h:54
bool global_tolerance
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accura...
Real max_deltat
Do not allow the adaptive time solver to select deltat > max_deltat.
Generic class from which first order UnsteadySolvers should subclass.
Data structure for holding completed parameter sensitivity calculations.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
Real upper_tolerance
This tolerance is the maximum relative error between an exact time integration and a single time step...
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with ...
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices) override=0
This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint s...
virtual void integrate_qoi_timestep() override=0
A method to integrate the system::QoI functionals.
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real min_deltat
Do not allow the adaptive time solver to select deltat < min_deltat.
virtual ~AdaptiveTimeSolver()
Destructor.
virtual void reinit() override
The reinitialization function.
virtual void retrieve_timestep() override
This method retrieves all the stored solutions at the current system.time.
virtual void advance_timestep() override
This method advances the solution to the next timestep, after a solve() has been performed.
virtual bool element_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
virtual void solve() override=0
This method solves for the solution at the next timestep.
std::unique_ptr< UnsteadySolver > core_time_solver
This object is used to take timesteps.
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver() override
An implicit linear solver to use for adjoint and sensitivity problems.
virtual void init() override
The initialization function.
Real max_growth
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old de...