Loading [MathJax]/extensions/tex2jax.js
Cantera  3.2.0a1
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
SystemJacobian.h
Go to the documentation of this file.
1//! @file SystemJacobian.h Declarations for class SystemJacobian
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
6#ifndef SYSTEMJACOBIAN_H
7#define SYSTEMJACOBIAN_H
8
10#include "cantera/base/global.h"
11
12namespace Cantera
13{
14
15/**
16 * Specifies the side of the system on which the preconditioner is applied. Not all
17 * methods are supported by all integrators.
18 */
20 NO_PRECONDITION, //! No preconditioning
21 LEFT_PRECONDITION, //! Left side preconditioning
22 RIGHT_PRECONDITION, //! Right side preconditioning
23 BOTH_PRECONDITION //! Left and right side preconditioning
24};
25
26//! Abstract base class representing Jacobian matrices and preconditioners used in
27//! nonlinear solvers.
28//!
29//! The matrix is initially populated with elements of the (steady-state) Jacobian
30//! using the setValue() method.
31//!
32//! If the object is being used for preconditioning, the preconditioner matrix
33//! @f$ M = I - \gamma J @f$ is computed and factorized by calling
34//! updatePreconditioner(). If the object is being used for solving a (pseudo-)transient
35//! problem, the transient Jacobian @f$ M = J - \alpha \Delta t^{-1} @f$ is computed and
36//! factorized by calling the updateTransient() method.
37//!
38//! Once the system is factorized, the solve() method can be used for preconditioning
39//! or solving Newton steps as appropriate.
40//!
41//! Implementations of this class provide different options for how to store the
42//! matrix elements and different algorithms for factorizing the matrix.
44{
45public:
47
48 virtual ~SystemJacobian() {}
49
50 //! Derived type, corresponding to names registered with SystemJacobianFactory
51 virtual const string type() const = 0;
52
53 //! Set a value at the specified row and column of the jacobian triplet vector
54 //! @param row row in the jacobian matrix
55 //! @param col column in the jacobian matrix
56 //! @param value value of the element at row and col
57 virtual void setValue(size_t row, size_t col, double value) {
58 throw NotImplementedError("SystemJacobian::setValue");
59 }
60
61 //! Adjust the state vector based on the preconditioner, e.g., Adaptive
62 //! preconditioning uses a strictly positive composition when preconditioning which
63 //! is handled by this function
64 //! @param state a vector containing the state to be updated
65 virtual void stateAdjustment(vector<double>& state) {
66 throw NotImplementedError("SystemJacobian::stateAdjustment");
67 }
68
69 //! Get preconditioner application side for CVODES
70 string preconditionerSide() const {
71 return m_precon_side;
72 }
73
74 //! For iterative solvers, set the side where the preconditioner is applied
75 virtual void setPreconditionerSide(const string& preconSide) {
76 m_precon_side = preconSide;
77 }
78
79 //! Transform Jacobian vector and write into
80 //! preconditioner, P = (I - gamma * J)
81 virtual void updatePreconditioner() {
82 throw NotImplementedError("SystemJacobian::updatePreconditioner");
83 }
84
85 //! Update the diagonal terms in the Jacobian by using the transient mask
86 //! @f$ \alpha @f$.
87 //! @param rdt Reciprocal of the time step [1/s]
88 //! @param mask Mask for transient terms: 1 if transient, 0 if algebraic.
89 virtual void updateTransient(double rdt, int* mask) {
90 throw NotImplementedError("SystemJacobian::updateTransient");
91 }
92
93 //! Solve a linear system using the system matrix *M*
94 //!
95 //! The matrix *M* can either be the preconditioner or the transient Jacobian
96 //! matrix, depending on whether updateTransient() or updatePreconditioner() was
97 //! called last.
98 //!
99 //! @param[in] stateSize length of the rhs and output vectors
100 //! @param[in] rhs_vector right hand side vector used in linear system
101 //! @param[out] output output vector for solution
102 virtual void solve(const size_t stateSize, double* rhs_vector, double* output) {
103 throw NotImplementedError("SystemJacobian::solve");
104 };
105
106 //! Perform preconditioner specific post-reactor setup operations such as factorize.
107 //! @deprecated To be removed after %Cantera 3.2. Use updatePreconditioner() and
108 //! factorize() instead.
109 virtual void setup() {
110 throw NotImplementedError("SystemJacobian::setup");
111 };
112
113 //! Reset parameters as needed
114 virtual void reset() {
115 throw NotImplementedError("SystemJacobian::reset");
116 };
117
118 //! Called during setup for any processes that need
119 //! to be completed prior to setup functions used in sundials.
120 //! @param networkSize the number of variables in the associated reactor network
121 virtual void initialize(size_t networkSize) {
122 throw NotImplementedError("SystemJacobian::initialize");
123 };
124
125 //! Used to provide system bandwidth for implementations that use banded matrix
126 //! storage. Ignored if not needed.
127 virtual void setBandwidth(size_t bw) {}
128
129 //! Print preconditioner contents
130 virtual void printPreconditioner() {
131 throw NotImplementedError("SystemJacobian::printPreconditioner");
132 };
133
134 //! Set gamma used in preconditioning
135 //! @param gamma used in M = I - gamma*J
136 virtual void setGamma(double gamma) {
137 m_gamma = gamma;
138 };
139
140 //! Get gamma used in preconditioning
141 virtual double gamma() {
142 return m_gamma;
143 };
144
145 //! Set the absolute tolerance in the solver outside of the network initialization
146 //! @param atol the specified tolerance
147 virtual void setAbsoluteTolerance(double atol) {
148 m_atol = atol;
149 }
150
151 //! Get latest status of linear solver. Zero indicates success. Meaning of non-zero
152 //! values is solver dependent.
153 virtual int info() const {
154 throw NotImplementedError("SystemJacobian::info");
155 }
156
157 //! Elapsed CPU time spent computing the Jacobian elements.
158 double elapsedTime() const {
159 return m_elapsed;
160 }
161
162 //! Increase the elapsed time by the specified time.
163 void updateElapsed(double evalTime) {
164 m_elapsed += evalTime;
165 }
166
167 //! Number of Jacobian evaluations.
168 int nEvals() const {
169 return m_nevals;
170 }
171 //! Increment the number of times the Jacobian has been evaluated
173 m_nevals++;
174 }
175
176 //! Number of times 'incrementAge' has been called since the last evaluation
177 int age() const {
178 return m_age;
179 }
180
181 //! Increment the Jacobian age.
183 m_age++;
184 }
185
186 //! Set the Jacobian age.
187 void setAge(int age) {
188 m_age = age;
189 }
190
191 //! Clear collected stats about number of Jacobian evaluations, CPU time spent on
192 //! Jacobian updates, and the number of Newton steps that have been taken since
193 //! updating the Jacobian.
194 void clearStats() {
195 m_elapsed = 0.0;
196 m_nevals = 0;
197 m_age = 100000;
198 }
199
200protected:
201 //! Factorize the system matrix. This method should be called at the end of
202 //! implementations of updatePreconditioner() and updateTransient().
203 virtual void factorize() {
204 throw NotImplementedError("SystemJacobian::factorize");
205 }
206
207 //! Dimension of the system
208 size_t m_dim;
209
210 //! gamma value used in M = I - gamma*J
211 double m_gamma = 1.0;
212
213 //! bool saying whether or not the system is initialized
214 bool m_init = false;
215
216 //! Absolute tolerance of the ODE solver
217 double m_atol = 0;
218
219 double m_elapsed = 0.0; //!< Elapsed CPU time taken to compute the Jacobian
220 int m_nevals = 0; //!< Number of Jacobian evaluations.
221
222 int m_age = 100000; //!< Age of the Jacobian (times incrementAge() has been called)
223
224 //! For iterative solvers, side of the system to apply the preconditioner on
225 string m_precon_side = "none";
226};
227
228//! @deprecated To be removed after Cantera 3.2. Renamed to SystemJacobian.
230{
232 warn_deprecated("PreconditionerBase::PreconditionerBase",
233 "To be removed after Cantera 3.2. Renamed to SystemJacobian.");
234 }
235};
236
237}
238#endif
An error indicating that an unimplemented function has been called.
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
virtual double gamma()
Get gamma used in preconditioning.
int nEvals() const
Number of Jacobian evaluations.
void setAge(int age)
Set the Jacobian age.
virtual void setGamma(double gamma)
Set gamma used in preconditioning.
size_t m_dim
Dimension of the system.
void incrementEvals()
Increment the number of times the Jacobian has been evaluated.
string preconditionerSide() const
Get preconditioner application side for CVODES.
virtual void setup()
Perform preconditioner specific post-reactor setup operations such as factorize.
double m_elapsed
Elapsed CPU time taken to compute the Jacobian.
double m_gamma
gamma value used in M = I - gamma*J
virtual int info() const
Get latest status of linear solver.
virtual void setPreconditionerSide(const string &preconSide)
For iterative solvers, set the side where the preconditioner is applied.
virtual void printPreconditioner()
Print preconditioner contents.
virtual void factorize()
Factorize the system matrix.
virtual void solve(const size_t stateSize, double *rhs_vector, double *output)
Solve a linear system using the system matrix M
virtual void reset()
Reset parameters as needed.
string m_precon_side
For iterative solvers, side of the system to apply the preconditioner on.
virtual void setBandwidth(size_t bw)
Used to provide system bandwidth for implementations that use banded matrix storage.
void incrementAge()
Increment the Jacobian age.
virtual void stateAdjustment(vector< double > &state)
Adjust the state vector based on the preconditioner, e.g., Adaptive preconditioning uses a strictly p...
void updateElapsed(double evalTime)
Increase the elapsed time by the specified time.
int age() const
Number of times 'incrementAge' has been called since the last evaluation.
virtual const string type() const =0
Derived type, corresponding to names registered with SystemJacobianFactory.
int m_nevals
Number of Jacobian evaluations.
bool m_init
bool saying whether or not the system is initialized
int m_age
Age of the Jacobian (times incrementAge() has been called)
void clearStats()
Clear collected stats about number of Jacobian evaluations, CPU time spent on Jacobian updates,...
virtual void initialize(size_t networkSize)
Called during setup for any processes that need to be completed prior to setup functions used in sund...
double m_atol
Absolute tolerance of the ODE solver.
virtual void setAbsoluteTolerance(double atol)
Set the absolute tolerance in the solver outside of the network initialization.
virtual void setValue(size_t row, size_t col, double value)
Set a value at the specified row and column of the jacobian triplet vector.
double elapsedTime() const
Elapsed CPU time spent computing the Jacobian elements.
virtual void updatePreconditioner()
Transform Jacobian vector and write into preconditioner, P = (I - gamma * J)
virtual void updateTransient(double rdt, int *mask)
Update the diagonal terms in the Jacobian by using the transient mask .
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
PreconditionerSide
Specifies the side of the system on which the preconditioner is applied.
@ LEFT_PRECONDITION
No preconditioning.
@ BOTH_PRECONDITION
Right side preconditioning.
@ RIGHT_PRECONDITION
Left side preconditioning.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Definition AnyMap.cpp:1997