https://mooseframework.inl.gov
TimeIntegrator.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 // MOOSE includes
13 #include "MooseObject.h"
14 #include "Restartable.h"
17 
18 class FEProblemBase;
19 class SystemBase;
21 
22 namespace libMesh
23 {
24 template <typename T>
25 class NumericVector;
26 } // namespace libMesh
27 
39 class TimeIntegrator : public MooseObject,
40  public Restartable,
43 {
44 public:
46 
48 
61  virtual void initialSetup() {}
62 
67  virtual void init();
68  virtual void preSolve() {}
69  virtual void preStep() {}
70 
74  virtual void solve();
75 
83  virtual void postSolve() {}
84 
88  virtual void postStep() {}
89 
90  virtual int order() = 0;
91 
100  virtual void computeTimeDerivatives() = 0;
101 
105  virtual void computeADTimeDerivatives(ADReal & ad_u_dot,
106  const dof_id_type & dof,
107  ADReal & ad_u_dot_dot) const = 0;
108 
112  virtual unsigned int getNumNonlinearIterations() const { return _n_nonlinear_iterations; }
113 
117  virtual unsigned int getNumLinearIterations() const { return _n_linear_iterations; }
118 
123  const Real & dt() const { return _dt; }
124 
128  virtual bool isExplicit() const { return false; }
129 
130  /*
131  * Returns whether the time integrator controls its own state. Explicit
132  * methods require extra care for determing when to store old solutions and
133  * stateful material properties.
134  */
135  virtual bool advancesProblemState() const { return false; }
136 
141  virtual unsigned int numStatesRequired() const { return 1; }
142 
146  virtual const bool & isLumped() const { return _is_lumped; }
147 
151  bool integratesVar(const unsigned int var_num) const;
152 
157 
161  virtual bool overridesSolve() const = 0;
162 
163 protected:
167  unsigned int getNumNonlinearIterationsLastSolve() const;
168 
172  unsigned int getNumLinearIterationsLastSolve() const;
173 
178  void copyVector(const NumericVector<Number> & from, NumericVector<Number> & to);
179 
184  virtual Real duDotDuCoeff() const { return 1; }
185 
189  void computeDuDotDu();
190 
193 
196 
200  std::vector<Real> & _du_dot_du;
201 
206  std::unique_ptr<NumericVector<Number>> & _solution_sub;
207  std::unique_ptr<NumericVector<Number>> & _solution_old_sub;
209 
211  int & _t_step;
212 
215 
218 
222  unsigned int _n_linear_iterations;
223 
226 
229 
232  std::vector<dof_id_type> & _local_indices;
233 
235  std::unordered_set<unsigned int> & _vars;
236 
238  std::unique_ptr<NumericVector<Number>> _from_subvector;
239 };
virtual void solve()
Solves the time step and sets the number of nonlinear and linear iterations.
virtual Real duDotDuCoeff() const
A class for creating restricted objects.
Definition: Restartable.h:28
virtual void preSolve()
bool & _var_restriction
Whether the user has requested that the time integrator be applied to a subset of variables...
void setNumIterationsLastSolve()
Record the linear and nonlinear iterations from a just finished solve.
virtual void postStep()
Callback to the TimeIntegrator called at the very end of time step.
std::unordered_set< unsigned int > & _vars
The variables that this time integrator integrates.
bool _is_lumped
Boolean flag that is set to true if lumped mass matrix is used.
FEProblemBase & _fe_problem
Reference to the problem.
Interface class for routines and member variables for time integrators relying on Newton&#39;s method...
virtual bool overridesSolve() const =0
SystemBase & _sys
Reference to the system this time integrator operates on.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
Base class for a system (of equations)
Definition: SystemBase.h:84
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
virtual unsigned int getNumNonlinearIterations() const
Gets the total number of nonlinear iterations over all stages of the time step.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
Real & _dt
The current time step size.
Real & _dt_old
The previous time step size.
virtual bool isExplicit() const
Returns whether the explicit solvers are used.
virtual void init()
Called only before the very first timestep (t_step = 0) Never called again (not even during recover/r...
Nonlinear system to be solved.
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
virtual void computeADTimeDerivatives(ADReal &ad_u_dot, const dof_id_type &dof, ADReal &ad_u_dot_dot) const =0
method for computing local automatic differentiation time derivatives
const Real & dt() const
Returns the time step size.
virtual unsigned int getNumLinearIterations() const
Gets the total number of linear iterations over all stages of the time step.
unsigned int _n_linear_iterations
Total number of linear iterations over all stages of the time step.
virtual int order()=0
std::vector< dof_id_type > & _local_indices
The local degree of freedom indices this time integrator is being applied to.
virtual void preStep()
std::unique_ptr< NumericVector< Number > > & _solution_old_sub
std::vector< Real > & _du_dot_du
Derivative of time derivative with respect to current solution: for the different variables...
virtual void initialSetup()
Called to setup datastructures.
Interface class for routines and member variables for time integrators relying on linear system assem...
virtual bool advancesProblemState() const
const NumericVector< Number > *const & _solution
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
unsigned int getNumLinearIterationsLastSolve() const
Gets the number of linear iterations in the most recent solve.
virtual const bool & isLumped() const
Returns whether mass matrix is lumped.
virtual void postSolve()
Callback to the TimeIntegrator called immediately after TimeIntegrator::solve() (so the name does mak...
const InputParameters & parameters() const
Get the parameters of the object.
const NumericVector< Number > & _solution_old
bool integratesVar(const unsigned int var_num) const
unsigned int _n_nonlinear_iterations
Total number of nonlinear iterations over all stages of the time step.
unsigned int getNumNonlinearIterationsLastSolve() const
Gets the number of nonlinear iterations in the most recent solve.
TimeIntegrator(const InputParameters &parameters)
int & _t_step
The current time step number.
std::unique_ptr< NumericVector< Number > > _from_subvector
A vector that is used for creating &#39;from&#39; subvectors in the copyVector() method.
static InputParameters validParams()
std::unique_ptr< NumericVector< Number > > & _solution_sub
uint8_t dof_id_type
virtual void computeTimeDerivatives()=0
Computes the time derivative and the Jacobian of the time derivative.
void computeDuDotDu()
Compute _du_dot_du.
void copyVector(const NumericVector< Number > &from, NumericVector< Number > &to)
Copy from one vector into another.
virtual unsigned int numStatesRequired() const
Return the number of states this requires in a linear system setting.