Cantera  3.1.0b1
Loading...
Searching...
No Matches
MultiJac.h
Go to the documentation of this file.
1//! @file MultiJac.h
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 CT_MULTIJAC_H
7#define CT_MULTIJAC_H
8
10#include "OneDim.h"
11
12namespace Cantera
13{
14
15/**
16 * Class MultiJac evaluates the Jacobian of a system of equations defined by a
17 * residual function supplied by an instance of class OneDim. The residual
18 * function may consist of several linked 1D domains, with different variables
19 * in each domain.
20 * @ingroup onedUtilsGroup
21 * @ingroup derivGroup
22 */
23class MultiJac : public BandMatrix
24{
25public:
26 //! Constructor.
27 //! @param r The nonlinear system for which to compute the Jacobian.
28 MultiJac(OneDim& r);
29
30 /**
31 * Evaluates the Jacobian at x0 using finite differences. The unperturbed residual
32 * function is resid0, which must be supplied on input. The parameter 'rdt' is the
33 * reciprocal of the time step. If zero, the steady-state Jacobian is evaluated,
34 * otherwise the transient Jacobian is evaluated.
35 *
36 * Each variable in the input vector `x0` is perturbed to compute the
37 * corresponding column of the Jacobian matrix. The Jacobian is computed by
38 * perturbing each variable, evaluating the residual function, and then
39 * estimating the partial derivatives numerically using finite differences.
40 *
41 * @param x0 Solution vector at which to evaluate the Jacobian
42 * @param resid0 Residual vector at x0
43 * @param rdt Reciprocal of the time step
44 */
45 void eval(double* x0, double* resid0, double rdt);
46
47 //! Elapsed CPU time spent computing the Jacobian.
48 double elapsedTime() const {
49 return m_elapsed;
50 }
51
52 //! Number of Jacobian evaluations.
53 int nEvals() const {
54 return m_nevals;
55 }
56
57 //! Number of times 'incrementAge' has been called since the last evaluation
58 int age() const {
59 return m_age;
60 }
61
62 //! Increment the Jacobian age.
63 void incrementAge() {
64 m_age++;
65 }
66
67 //! Update the transient terms in the Jacobian by using the transient mask.
68 void updateTransient(double rdt, integer* mask);
69
70 //! Set the Jacobian age.
71 void setAge(int age) {
72 m_age = age;
73 }
74
75 //! Return the transient mask.
76 vector<int>& transientMask() {
77 return m_mask;
78 }
79
80 //! @deprecated To be removed after Cantera 3.1.
81 void incrementDiagonal(int j, double d);
82
83protected:
84 //! Residual evaluator for this Jacobian
85 /*!
86 * This is a pointer to the residual evaluator. This object isn't owned by
87 * this Jacobian object.
88 */
90
91 vector<double> m_r1; //!< Perturbed residual vector
92 double m_rtol = 1e-5; //!< Relative tolerance for perturbing solution components
93
94 //! Absolute tolerance for perturbing solution components
95 double m_atol = sqrt(std::numeric_limits<double>::epsilon());
96
97 double m_elapsed = 0.0; //!< Elapsed CPU time taken to compute the Jacobian
98 vector<double> m_ssdiag; //!< Diagonal of the steady-state Jacobian
99
100 //! Transient mask for transient terms, 1 if transient, 0 if steady-state
101 vector<int> m_mask;
102
103 int m_nevals = 0; //!< Number of Jacobian evaluations.
104 int m_age = 100000; //!< Age of the Jacobian (times incrementAge() has been called)
105};
106}
107
108#endif
Declarations for the class BandMatrix which is a child class of GeneralMatrix for banded matrices han...
A class for banded matrices, involving matrix inversion processes.
Definition BandMatrix.h:37
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Definition MultiJac.h:24
double m_rtol
Relative tolerance for perturbing solution components.
Definition MultiJac.h:92
int nEvals() const
Number of Jacobian evaluations.
Definition MultiJac.h:53
void setAge(int age)
Set the Jacobian age.
Definition MultiJac.h:71
vector< int > & transientMask()
Return the transient mask.
Definition MultiJac.h:76
double m_elapsed
Elapsed CPU time taken to compute the Jacobian.
Definition MultiJac.h:97
void incrementDiagonal(int j, double d)
Definition MultiJac.cpp:28
vector< int > m_mask
Transient mask for transient terms, 1 if transient, 0 if steady-state.
Definition MultiJac.h:101
void eval(double *x0, double *resid0, double rdt)
Evaluates the Jacobian at x0 using finite differences.
Definition MultiJac.cpp:35
void updateTransient(double rdt, integer *mask)
Update the transient terms in the Jacobian by using the transient mask.
Definition MultiJac.cpp:21
void incrementAge()
Increment the Jacobian age.
Definition MultiJac.h:63
int age() const
Number of times 'incrementAge' has been called since the last evaluation.
Definition MultiJac.h:58
vector< double > m_ssdiag
Diagonal of the steady-state Jacobian.
Definition MultiJac.h:98
int m_nevals
Number of Jacobian evaluations.
Definition MultiJac.h:103
int m_age
Age of the Jacobian (times incrementAge() has been called)
Definition MultiJac.h:104
OneDim * m_resid
Residual evaluator for this Jacobian.
Definition MultiJac.h:89
vector< double > m_r1
Perturbed residual vector.
Definition MultiJac.h:91
double m_atol
Absolute tolerance for perturbing solution components.
Definition MultiJac.h:95
double elapsedTime() const
Elapsed CPU time spent computing the Jacobian.
Definition MultiJac.h:48
Container class for multiple-domain 1D problems.
Definition OneDim.h:27
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595