https://mooseframework.inl.gov
ExplicitTimeIntegrator.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 #include "TimeIntegrator.h"
13 #include "MeshChangedInterface.h"
14 
15 #include "libmesh/linear_solver.h"
16 #include "libmesh/preconditioner.h"
17 #include "libmesh/numeric_vector.h"
18 
19 #include "LumpedPreconditioner.h"
20 
26 {
27 public:
29 
31 
32  virtual void initialSetup() override;
33  virtual void init() override;
34  virtual void preSolve() override;
35  virtual void meshChanged() override;
36  virtual bool isExplicit() const override { return true; }
37 
38 protected:
39  enum SolveType
40  {
44  };
45 
51  virtual bool performExplicitSolve(SparseMatrix<Number> & mass_matrix);
52 
58  bool solveLinearSystem(SparseMatrix<Number> & mass_matrix);
59 
64 
68  virtual TagID massMatrixTagID() const { return _Ke_time_tag; }
69 
72 
74  NumericVector<Real> * _explicit_residual;
75 
77  NumericVector<Real> * _solution_update;
78 
80  NumericVector<Real> * _mass_matrix_diag_inverted;
81 
83  NumericVector<Real> * _ones;
84 
87 
89  std::unique_ptr<libMesh::LinearSolver<Number>> _linear_solver;
90 
92  std::unique_ptr<LumpedPreconditioner> _preconditioner;
93 
96 };
bool checkLinearConvergence()
Check for the linear solver convergence.
NumericVector< Real > * _solution_update
Solution vector for the linear solve.
Real _current_time
Save off current time to reset it back and forth.
unsigned int TagID
Definition: MooseTypes.h:210
NumericVector< Real > * _ones
Vector of 1&#39;s to help with creating the lumped mass matrix.
virtual void initialSetup() override
Called to setup datastructures.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:127
TagID _Ke_time_tag
For computing the mass matrix.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
MooseEnum _solve_type
Solve type for how mass matrix is handled.
static InputParameters validParams()
Interface for notifications that the mesh has changed.
std::unique_ptr< libMesh::LinearSolver< Number > > _linear_solver
For solving with the consistent matrix.
NumericVector< Real > * _explicit_residual
Residual used for the RHS.
virtual void preSolve() override
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
virtual bool performExplicitSolve(SparseMatrix< Number > &mass_matrix)
Solves a linear system using the chosen solve type.
std::unique_ptr< LumpedPreconditioner > _preconditioner
For solving with lumped preconditioning.
NumericVector< Real > * _mass_matrix_diag_inverted
Diagonal of the lumped mass matrix (and its inversion)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
virtual void init() override
Called only before the very first timestep (t_step = 0) Never called again (not even during recover/r...
virtual TagID massMatrixTagID() const
virtual void meshChanged() override
Called on this object when the mesh changes.
virtual bool isExplicit() const override
Returns whether the explicit solvers are used.
bool solveLinearSystem(SparseMatrix< Number > &mass_matrix)
Solves a linear system.
Base class for explicit time integrators that are implemented without using a nonlinear solver...
ExplicitTimeIntegrator(const InputParameters &parameters)