15 #include "libmesh/transient_system.h" 16 #include "libmesh/linear_implicit_system.h" 17 #include "libmesh/linear_solver.h" 41 virtual void solve()
override;
61 const std::set<TagID> & vector_tags_to_close)
override;
80 const std::set<TagID> & matrix_tags,
81 const bool compute_gradients =
true);
99 std::vector<dof_id_type> & n_nz,
100 std::vector<dof_id_type> & n_oz)
override;
108 virtual const System &
system()
const override {
return _sys; }
153 const std::set<TagID> & matrix_tags,
154 const bool compute_gradients =
true);
virtual void stopSolve(const ExecFlagType &exec_flag, const std::set< TagID > &vector_tags_to_close) override
Quit the current solve as soon as possible.
Base class for finite volume kernels that contribute to a linear systems.
void computeGradients()
Compute the Green-Gauss gradients.
unsigned int _n_linear_iters
Number of linear iterations.
TagID _system_matrix_non_time_tag
Tag for non-time contribution to the system matrix.
TagID rightHandSideVectorTag() const
virtual bool containsTimeKernel() override
If the system has a kernel that corresponds to a time derivative.
libMesh::LinearImplicitSystem & linearImplicitSystem()
Return a reference to the stored linear implicit system.
std::vector< std::unique_ptr< NumericVector< Number > > > & newGradientContainer()
Return a reference to the new (temporary) gradient container vectors.
NumericVector< Number > & getRightHandSideNonTimeVector()
Return a numeric vector that is associated with the nontime tag.
virtual TagID systemMatrixTag() const override
Return the Matrix Tag ID for System.
const NumericVector< Number > & getRightHandSideVector() const
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.
virtual void compute(ExecFlagType type) override
Compute time derivatives, auxiliary variables, etc.
TagID rightHandSideNonTimeVectorTag() const
NumericVector< Number > * rhs
std::set< TagID > _vector_tags
Vector tags to temporarily store all tags associated with the current system.
TagID rightHandSideTimeVectorTag() const
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
unsigned int nLinearIterations() const
Return the number of linear iterations.
unsigned int _current_l_its
The linear iterations needed for convergence.
SparseMatrix< Number > & getSystemMatrix()
Fetching the system matrix from the libmesh system.
Real _initial_linear_residual
The initial linear residual.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
const SparseMatrix< Number > & getSystemMatrix() const
std::vector< std::unique_ptr< NumericVector< Number > > > _new_gradient
Vectors to store the new gradients during the computation.
System & _sys
Base class reference to the libmesh system.
virtual void initialSetup() override
Setup Functions.
virtual const std::string & name() const
libMesh::LinearImplicitSystem & _linear_implicit_system
Base class reference to the linear implicit system in libmesh.
NumericVector< Number > * _rhs_time
right hand side vector for time contributions
Real _final_linear_residual
The final linear residual.
std::set< TagID > _matrix_tags
Matrix tags to temporarily store all tags associated with the current system.
std::vector< NumericVector< Number > * > _solution_state
The current states of the solution (0 = current, 1 = old, etc)
TagID _rhs_tag
Used for the right hand side vector from PETSc.
virtual void jacobianSetup() override
TagID _rhs_non_time_tag
Tag for non-time contribution rhs.
virtual void solve() override
Solve the system (using libMesh magic)
virtual const System & system() const override
Interface for objects interacting with the PerfGraph.
TagID _rhs_time_tag
Tag for time contribution rhs.
TagID _system_matrix_tag
Tag for every contribution to system matrix.
bool _converged
If the solve on the linear system converged.
NumericVector< Number > * _rhs_non_time
right hand side vector for non-time contributions
NumericVector< Number > & getRightHandSideVector()
Fetching the right hand side vector from the libmesh system.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Class for containing MooseEnum item information.
NumericVector< Number > & getRightHandSideTimeVector()
Return a numeric vector that is associated with the time tag.
SparseMatrix< Number > * matrix
virtual void residualSetup() override
virtual std::vector< std::string > timeKernelVariableNames() override
Returns the names of the variables that have time derivative kernels in the system.
Linear system to be solved.
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.
LinearSystem(FEProblemBase &problem, const std::string &name)
virtual System & system() override
Get the reference to the libMesh system.
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.
virtual bool converged() override
At the moment, this is only used for the multi-system fixed point iteration.