https://mooseframework.inl.gov
LinearSystem.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 "SolverSystem.h"
13 #include "PerfGraphInterface.h"
14 
15 #include "libmesh/transient_system.h"
16 #include "libmesh/linear_implicit_system.h"
17 #include "libmesh/linear_solver.h"
18 
19 class LinearFVKernel;
20 
21 // libMesh forward declarations
22 namespace libMesh
23 {
24 template <typename T>
25 class NumericVector;
26 template <typename T>
27 class SparseMatrix;
28 template <typename T>
29 class DiagonalMatrix;
30 } // namespace libMesh
31 
36 {
37 public:
38  LinearSystem(FEProblemBase & problem, const std::string & name);
39  virtual ~LinearSystem();
40 
41  virtual void solve() override;
42 
48  virtual bool converged() override { return _converged; }
49 
50  virtual void initialSetup() override;
51 
52  // Overriding these to make sure the linear systems don't do anything during
53  // residual/jacobian setup
54  virtual void residualSetup() override {}
55  virtual void jacobianSetup() override {}
56 
60  virtual void stopSolve(const ExecFlagType & exec_flag,
61  const std::set<TagID> & vector_tags_to_close) override;
62 
68  virtual bool containsTimeKernel() override;
69  virtual std::vector<std::string> timeKernelVariableNames() override { return {}; }
70 
79  void computeLinearSystemTags(const std::set<TagID> & vector_tags,
80  const std::set<TagID> & matrix_tags,
81  const bool compute_gradients = true);
82 
87 
92 
97 
98  virtual void augmentSparsity(SparsityPattern::Graph & sparsity,
99  std::vector<dof_id_type> & n_nz,
100  std::vector<dof_id_type> & n_oz) override;
101 
105  unsigned int nLinearIterations() const { return _n_linear_iters; }
106 
107  virtual System & system() override { return _sys; }
108  virtual const System & system() const override { return _sys; }
109 
115  virtual TagID systemMatrixTag() const override { return _system_matrix_tag; }
117 
121  {
123  }
124 
128 
132  void computeGradients();
133 
137  std::vector<std::unique_ptr<NumericVector<Number>>> & newGradientContainer()
138  {
139  return _new_gradient;
140  }
141 
142  virtual void compute(ExecFlagType type) override;
143 
144 protected:
152  void computeLinearSystemInternal(const std::set<TagID> & vector_tags,
153  const std::set<TagID> & matrix_tags,
154  const bool compute_gradients = true);
155 
157  System & _sys;
158 
160  unsigned int _current_l_its;
161 
163  std::set<TagID> _vector_tags;
164 
166  std::set<TagID> _matrix_tags;
167 
170 
173 
176 
179 
182 
185 
188 
190  unsigned int _n_linear_iters;
191 
194 
197 
200 
203 
208  std::vector<std::unique_ptr<NumericVector<Number>>> _new_gradient;
209 
210 private:
212  std::vector<NumericVector<Number> *> _solution_state;
213 };
virtual void stopSolve(const ExecFlagType &exec_flag, const std::set< TagID > &vector_tags_to_close) override
Quit the current solve as soon as possible.
Definition: LinearSystem.C:301
Base class for finite volume kernels that contribute to a linear systems.
void computeGradients()
Compute the Green-Gauss gradients.
Definition: LinearSystem.C:142
unsigned int _n_linear_iters
Number of linear iterations.
Definition: LinearSystem.h:190
unsigned int TagID
Definition: MooseTypes.h:210
TagID _system_matrix_non_time_tag
Tag for non-time contribution to the system matrix.
Definition: LinearSystem.h:184
TagID rightHandSideVectorTag() const
Definition: LinearSystem.h:114
virtual bool containsTimeKernel() override
If the system has a kernel that corresponds to a time derivative.
Definition: LinearSystem.C:310
libMesh::LinearImplicitSystem & linearImplicitSystem()
Return a reference to the stored linear implicit system.
Definition: LinearSystem.h:86
std::vector< std::unique_ptr< NumericVector< Number > > > & newGradientContainer()
Return a reference to the new (temporary) gradient container vectors.
Definition: LinearSystem.h:137
NumericVector< Number > & getRightHandSideNonTimeVector()
Return a numeric vector that is associated with the nontime tag.
Definition: LinearSystem.C:260
virtual TagID systemMatrixTag() const override
Return the Matrix Tag ID for System.
Definition: LinearSystem.h:115
const NumericVector< Number > & getRightHandSideVector() const
Definition: LinearSystem.h:120
void computeLinearSystemInternal(const std::set< TagID > &vector_tags, const std::set< TagID > &matrix_tags, const bool compute_gradients=true)
Compute the right hand side and system matrix for given tags.
Definition: LinearSystem.C:185
virtual void compute(ExecFlagType type) override
Compute time derivatives, auxiliary variables, etc.
Definition: LinearSystem.C:332
TagID rightHandSideNonTimeVectorTag() const
Definition: LinearSystem.h:113
NumericVector< Number > * rhs
std::set< TagID > _vector_tags
Vector tags to temporarily store all tags associated with the current system.
Definition: LinearSystem.h:163
TagID rightHandSideTimeVectorTag() const
Definition: LinearSystem.h:112
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
unsigned int nLinearIterations() const
Return the number of linear iterations.
Definition: LinearSystem.h:105
unsigned int _current_l_its
The linear iterations needed for convergence.
Definition: LinearSystem.h:160
SparseMatrix< Number > & getSystemMatrix()
Fetching the system matrix from the libmesh system.
Definition: LinearSystem.h:126
Real _initial_linear_residual
The initial linear residual.
Definition: LinearSystem.h:193
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
const SparseMatrix< Number > & getSystemMatrix() const
Definition: LinearSystem.h:127
std::vector< std::unique_ptr< NumericVector< Number > > > _new_gradient
Vectors to store the new gradients during the computation.
Definition: LinearSystem.h:208
System & _sys
Base class reference to the libmesh system.
Definition: LinearSystem.h:157
virtual void initialSetup() override
Setup Functions.
Definition: LinearSystem.C:103
virtual const std::string & name() const
Definition: SystemBase.C:1330
libMesh::LinearImplicitSystem & _linear_implicit_system
Base class reference to the linear implicit system in libmesh.
Definition: LinearSystem.h:202
NumericVector< Number > * _rhs_time
right hand side vector for time contributions
Definition: LinearSystem.h:172
Real _final_linear_residual
The final linear residual.
Definition: LinearSystem.h:196
std::set< TagID > _matrix_tags
Matrix tags to temporarily store all tags associated with the current system.
Definition: LinearSystem.h:166
std::vector< NumericVector< Number > * > _solution_state
The current states of the solution (0 = current, 1 = old, etc)
Definition: LinearSystem.h:212
TagID _rhs_tag
Used for the right hand side vector from PETSc.
Definition: LinearSystem.h:181
virtual void jacobianSetup() override
Definition: LinearSystem.h:55
TagID _rhs_non_time_tag
Tag for non-time contribution rhs.
Definition: LinearSystem.h:175
virtual void solve() override
Solve the system (using libMesh magic)
Definition: LinearSystem.C:274
virtual const System & system() const override
Definition: LinearSystem.h:108
Interface for objects interacting with the PerfGraph.
TagID _rhs_time_tag
Tag for time contribution rhs.
Definition: LinearSystem.h:169
TagID _system_matrix_tag
Tag for every contribution to system matrix.
Definition: LinearSystem.h:187
bool _converged
If the solve on the linear system converged.
Definition: LinearSystem.h:199
NumericVector< Number > * _rhs_non_time
right hand side vector for non-time contributions
Definition: LinearSystem.h:178
NumericVector< Number > & getRightHandSideVector()
Fetching the right hand side vector from the libmesh system.
Definition: LinearSystem.h:119
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Class for containing MooseEnum item information.
Definition: MooseEnumItem.h:18
NumericVector< Number > & getRightHandSideTimeVector()
Return a numeric vector that is associated with the time tag.
Definition: LinearSystem.C:254
SparseMatrix< Number > * matrix
virtual void residualSetup() override
Definition: LinearSystem.h:54
virtual std::vector< std::string > timeKernelVariableNames() override
Returns the names of the variables that have time derivative kernels in the system.
Definition: LinearSystem.h:69
Linear system to be solved.
Definition: LinearSystem.h:35
virtual void augmentSparsity(SparsityPattern::Graph &sparsity, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz) override
Will modify the sparsity pattern to add logical geometric connections.
Definition: LinearSystem.C:266
virtual ~LinearSystem()
LinearSystem(FEProblemBase &problem, const std::string &name)
Definition: LinearSystem.C:75
virtual System & system() override
Get the reference to the libMesh system.
Definition: LinearSystem.h:107
void computeLinearSystemTags(const std::set< TagID > &vector_tags, const std::set< TagID > &matrix_tags, const bool compute_gradients=true)
Compute the right hand side and the system matrix of the system for given tags.
Definition: LinearSystem.C:116
virtual bool converged() override
At the moment, this is only used for the multi-system fixed point iteration.
Definition: LinearSystem.h:48