www.mooseframework.org
ExplicitTimeIntegrator.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
21 
22 template <>
24 
30 {
31 public:
33 
35 
36  virtual void initialSetup() override;
37  virtual void init() override;
38  virtual void preSolve() override;
39  virtual void meshChanged() override;
40 
41 protected:
42  enum SolveType
43  {
47  };
48 
54  bool performExplicitSolve(SparseMatrix<Number> & mass_matrix);
55 
61  bool solveLinearSystem(SparseMatrix<Number> & mass_matrix);
62 
67 
70 
72  NumericVector<Real> & _explicit_residual;
73 
75  NumericVector<Real> & _solution_update;
76 
78  NumericVector<Real> & _mass_matrix_diag;
79 
81  NumericVector<Real> * _ones;
82 
85 
87  std::unique_ptr<LinearSolver<Number>> _linear_solver;
88 
90  std::unique_ptr<LumpedPreconditioner> _preconditioner;
91 
94 };
95 
99 class LumpedPreconditioner : public Preconditioner<Real>
100 {
101 public:
102  static InputParameters validParams();
103 
104  LumpedPreconditioner(const NumericVector<Real> & diag_inverse)
105  : Preconditioner(diag_inverse.comm()), _diag_inverse(diag_inverse)
106  {
107  }
108 
109  virtual void init() override
110  {
111  // No more initialization needed here
112  _is_initialized = true;
113  }
114 
115  virtual void apply(const NumericVector<Real> & x, NumericVector<Real> & y) override
116  {
117  y.pointwise_mult(_diag_inverse, x);
118  }
119 
120 protected:
122  const NumericVector<Real> & _diag_inverse;
123 };
validParams< ExplicitTimeIntegrator >
InputParameters validParams< ExplicitTimeIntegrator >()
ExplicitTimeIntegrator::validParams
static InputParameters validParams()
Definition: ExplicitTimeIntegrator.C:22
comm
MPI_Comm comm
Definition: PetscDMMoose.C:1505
ExplicitTimeIntegrator::_Ke_time_tag
TagID _Ke_time_tag
For computing the mass matrix.
Definition: ExplicitTimeIntegrator.h:84
MeshChangedInterface.h
ExplicitTimeIntegrator::init
virtual void init() override
Called only before the very first timestep (t_step = 0) Never called again (not even during recover/r...
Definition: ExplicitTimeIntegrator.C:64
ExplicitTimeIntegrator::SolveType
SolveType
Definition: ExplicitTimeIntegrator.h:42
MooseEnum
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
ExplicitTimeIntegrator::_solve_type
MooseEnum _solve_type
Solve type for how mass matrix is handled.
Definition: ExplicitTimeIntegrator.h:69
TagID
unsigned int TagID
Definition: MooseTypes.h:197
TimeIntegrator
Base class for time integrators.
Definition: TimeIntegrator.h:43
MooseObject::parameters
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:76
ExplicitTimeIntegrator::CONSISTENT
Definition: ExplicitTimeIntegrator.h:44
ExplicitTimeIntegrator::_mass_matrix_diag
NumericVector< Real > & _mass_matrix_diag
Diagonal of the lumped mass matrix (and its inversion)
Definition: ExplicitTimeIntegrator.h:78
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
x
static PetscErrorCode Vec x
Definition: PetscDMMoose.C:1263
ExplicitTimeIntegrator::performExplicitSolve
bool performExplicitSolve(SparseMatrix< Number > &mass_matrix)
Solves a linear system using the chosen solve type.
Definition: ExplicitTimeIntegrator.C:95
ExplicitTimeIntegrator
Base class for explicit time integrators that are implemented without using a nonlinear solver.
Definition: ExplicitTimeIntegrator.h:29
LumpedPreconditioner::validParams
static InputParameters validParams()
ExplicitTimeIntegrator::_solution_update
NumericVector< Real > & _solution_update
Solution vector for the linear solve.
Definition: ExplicitTimeIntegrator.h:75
ExplicitTimeIntegrator::preSolve
virtual void preSolve() override
Definition: ExplicitTimeIntegrator.C:69
LumpedPreconditioner
Helper class to apply lumped mass matrix preconditioner.
Definition: ExplicitTimeIntegrator.h:99
ExplicitTimeIntegrator::_linear_solver
std::unique_ptr< LinearSolver< Number > > _linear_solver
For solving with the consistent matrix.
Definition: ExplicitTimeIntegrator.h:87
ExplicitTimeIntegrator::_explicit_residual
NumericVector< Real > & _explicit_residual
Residual used for the RHS.
Definition: ExplicitTimeIntegrator.h:72
LumpedPreconditioner::apply
virtual void apply(const NumericVector< Real > &x, NumericVector< Real > &y) override
Definition: ExplicitTimeIntegrator.h:115
MeshChangedInterface
Interface for notifications that the mesh has changed.
Definition: MeshChangedInterface.h:28
LumpedPreconditioner::LumpedPreconditioner
LumpedPreconditioner(const NumericVector< Real > &diag_inverse)
Definition: ExplicitTimeIntegrator.h:104
LumpedPreconditioner::_diag_inverse
const NumericVector< Real > & _diag_inverse
The inverse of the diagonal of the lumped matrix.
Definition: ExplicitTimeIntegrator.h:122
TimeIntegrator.h
ExplicitTimeIntegrator::_preconditioner
std::unique_ptr< LumpedPreconditioner > _preconditioner
For solving with lumped preconditioning.
Definition: ExplicitTimeIntegrator.h:90
ExplicitTimeIntegrator::LUMP_PRECONDITIONED
Definition: ExplicitTimeIntegrator.h:46
LumpedPreconditioner::init
virtual void init() override
Definition: ExplicitTimeIntegrator.h:109
ExplicitTimeIntegrator::solveLinearSystem
bool solveLinearSystem(SparseMatrix< Number > &mass_matrix)
Solves a linear system.
Definition: ExplicitTimeIntegrator.C:146
ExplicitTimeIntegrator::_current_time
Real _current_time
Save off current time to reset it back and forth.
Definition: ExplicitTimeIntegrator.h:93
ExplicitTimeIntegrator::ExplicitTimeIntegrator
ExplicitTimeIntegrator(const InputParameters &parameters)
Definition: ExplicitTimeIntegrator.C:39
ExplicitTimeIntegrator::meshChanged
virtual void meshChanged() override
Called on this object when the mesh changes.
Definition: ExplicitTimeIntegrator.C:74
ExplicitTimeIntegrator::_ones
NumericVector< Real > * _ones
Vector of 1's to help with creating the lumped mass matrix.
Definition: ExplicitTimeIntegrator.h:81
ExplicitTimeIntegrator::initialSetup
virtual void initialSetup() override
Called to setup datastructures.
Definition: ExplicitTimeIntegrator.C:58
ExplicitTimeIntegrator::checkLinearConvergence
bool checkLinearConvergence()
Check for the linear solver convergence.
Definition: ExplicitTimeIntegrator.C:165
ExplicitTimeIntegrator::LUMPED
Definition: ExplicitTimeIntegrator.h:45