Cantera  3.0.0
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 MultiJac(OneDim& r);
27
28 /**
29 * Evaluate the Jacobian at x0. The unperturbed residual function is resid0,
30 * which must be supplied on input. The third parameter 'rdt' is the
31 * reciprocal of the time step. If zero, the steady-state Jacobian is
32 * evaluated.
33 */
34 void eval(double* x0, double* resid0, double rdt);
35
36 //! Elapsed CPU time spent computing the Jacobian.
37 double elapsedTime() const {
38 return m_elapsed;
39 }
40
41 //! Number of Jacobian evaluations.
42 int nEvals() const {
43 return m_nevals;
44 }
45
46 //! Number of times 'incrementAge' has been called since the last evaluation
47 int age() const {
48 return m_age;
49 }
50
51 //! Increment the Jacobian age.
52 void incrementAge() {
53 m_age++;
54 }
55
56 void updateTransient(double rdt, integer* mask);
57
58 //! Set the Jacobian age.
59 void setAge(int age) {
60 m_age = age;
61 }
62
63 vector<int>& transientMask() {
64 return m_mask;
65 }
66
67 void incrementDiagonal(int j, double d);
68
69protected:
70 //! Residual evaluator for this Jacobian
71 /*!
72 * This is a pointer to the residual evaluator. This object isn't owned by
73 * this Jacobian object.
74 */
76
77 vector<double> m_r1;
78 double m_rtol = 1e-5;
79 double m_atol = sqrt(std::numeric_limits<double>::epsilon());
80 double m_elapsed = 0.0;
81 vector<double> m_ssdiag;
82 vector<int> m_mask;
83 int m_nevals = 0;
84 int m_age = 100000;
85 size_t m_size;
86 size_t m_points;
87};
88}
89
90#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
int nEvals() const
Number of Jacobian evaluations.
Definition MultiJac.h:42
void setAge(int age)
Set the Jacobian age.
Definition MultiJac.h:59
void eval(double *x0, double *resid0, double rdt)
Evaluate the Jacobian at x0.
Definition MultiJac.cpp:36
void incrementAge()
Increment the Jacobian age.
Definition MultiJac.h:52
int age() const
Number of times 'incrementAge' has been called since the last evaluation.
Definition MultiJac.h:47
OneDim * m_resid
Residual evaluator for this Jacobian.
Definition MultiJac.h:75
double elapsedTime() const
Elapsed CPU time spent computing the Jacobian.
Definition MultiJac.h:37
Container class for multiple-domain 1D problems.
Definition OneDim.h:27
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564