Cantera  3.3.0a1
Loading...
Searching...
No Matches
AdaptivePreconditioner.h
Go to the documentation of this file.
1/**
2 * @file AdaptivePreconditioner.h Declarations for the class
3 * AdaptivePreconditioner which is a child class of SystemJacobian
4 * for preconditioners used by sundials
5 */
6
7// This file is part of Cantera. See License.txt in the top-level directory or
8// at https://cantera.org/license.txt for license and copyright information.
9
10#ifndef ADAPTIVEPRECONDITIONER_H
11#define ADAPTIVEPRECONDITIONER_H
12
14
15namespace Cantera
16{
17
18//! AdaptivePreconditioner a preconditioner designed for use with large
19//! mechanisms that leverages sparse solvers. It does this by pruning
20//! the preconditioner by a threshold value. It also neglects pressure
21//! dependence and third body contributions in its formation and has a
22//! finite difference approximation for temperature.
24{
25public:
27
28 void initialize(size_t networkSize) override;
29 const string type() const override { return "Adaptive"; }
30 void factorize() override;
31 void solve(const size_t stateSize, double* rhs_vector, double* output) override;
32 void stateAdjustment(vector<double>& state) override;
33
34 int info() const override {
35 return static_cast<int>(m_solver.info());
36 }
37
38 //! Prune preconditioner elements
40
41 //! Get the threshold value for setting elements
42 double threshold() { return m_threshold; }
43
44 //! Get ILUT fill factor
45 double ilutFillFactor() { return m_fill_factor; }
46
47 //! Get ILUT drop tolerance
48 double ilutDropTol() { return m_drop_tol; }
49
50 //! Set the threshold value to compare elements against
51 //! @param threshold double value used in setting by threshold
52 void setThreshold(double threshold) {
54 m_prune_precon = (threshold <= 0) ? false : true;
55 }
56
57 //! Set drop tolerance for ILUT
58 //! @param droptol double value used in setting solver drop tolerance
59 void setIlutDropTol(double droptol) {
60 m_drop_tol = droptol;
61 m_solver.setDroptol(droptol);
62 }
63
64 //! Set the fill factor for ILUT solver
65 //! @param fillFactor fill in factor for ILUT solver
66 void setIlutFillFactor(int fillFactor) {
67 m_fill_factor = fillFactor;
68 m_solver.setFillfactor(fillFactor);
69 }
70
71protected:
72 //! ILUT fill factor
73 double m_fill_factor = 0;
74
75 //! ILUT drop tolerance
76 double m_drop_tol = 0;
77
78 //! Solver used in solving the linear system
79 Eigen::IncompleteLUT<double> m_solver;
80
81 //! Minimum value a non-diagonal element must be to be included in
82 //! the preconditioner
83 double m_threshold = 0.0;
84
85 //! Bool set whether to prune the matrix or not
86 double m_prune_precon = true;
87};
88
89}
90
91#endif
AdaptivePreconditioner a preconditioner designed for use with large mechanisms that leverages sparse ...
double threshold()
Get the threshold value for setting elements.
Eigen::IncompleteLUT< double > m_solver
Solver used in solving the linear system.
double ilutDropTol()
Get ILUT drop tolerance.
void setIlutDropTol(double droptol)
Set drop tolerance for ILUT.
double m_drop_tol
ILUT drop tolerance.
void solve(const size_t stateSize, double *rhs_vector, double *output) override
Solve a linear system using the system matrix M
double m_prune_precon
Bool set whether to prune the matrix or not.
double m_threshold
Minimum value a non-diagonal element must be to be included in the preconditioner.
void factorize() override
Factorize the system matrix.
void initialize(size_t networkSize) override
Called during setup for any processes that need to be completed prior to setup functions used in sund...
void setIlutFillFactor(int fillFactor)
Set the fill factor for ILUT solver.
void stateAdjustment(vector< double > &state) override
Adjust the state vector based on the preconditioner, e.g., Adaptive preconditioning uses a strictly p...
const string type() const override
Derived type, corresponding to names registered with SystemJacobianFactory.
void prunePreconditioner()
Prune preconditioner elements.
void setThreshold(double threshold)
Set the threshold value to compare elements against.
int info() const override
Get latest status of linear solver.
double ilutFillFactor()
Get ILUT fill factor.
System Jacobians that use Eigen sparse matrices for storage.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595