Loading [MathJax]/extensions/tex2jax.js
Cantera  3.2.0a1
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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
11#include "OneDim.h"
12
13namespace Cantera
14{
15
16/**
17 * Class MultiJac evaluates the Jacobian of a system of equations defined by a
18 * residual function supplied by an instance of class OneDim. The residual
19 * function may consist of several linked 1D domains, with different variables
20 * in each domain.
21 * @ingroup onedUtilsGroup
22 * @ingroup derivGroup
23 */
25{
26public:
27 //! Constructor.
28 //! @param r The nonlinear system for which to compute the Jacobian.
29 //! @deprecated To be removed after %Cantera 3.2. Use default constructor instead.
30 MultiJac(OneDim& r);
31
32 MultiJac() = default;
33
34 /**
35 * Evaluates the Jacobian at x0 using finite differences. The unperturbed residual
36 * function is resid0, which must be supplied on input. The parameter 'rdt' is the
37 * reciprocal of the time step. If zero, the steady-state Jacobian is evaluated,
38 * otherwise the transient Jacobian is evaluated.
39 *
40 * Each variable in the input vector `x0` is perturbed to compute the
41 * corresponding column of the Jacobian matrix. The Jacobian is computed by
42 * perturbing each variable, evaluating the residual function, and then
43 * estimating the partial derivatives numerically using finite differences.
44 *
45 * @param x0 Solution vector at which to evaluate the Jacobian
46 * @param resid0 Residual vector at x0
47 * @param rdt Reciprocal of the time step
48 * @deprecated To be removed after %Cantera 3.2. Jacobian evaluation moved to
49 * OneDim::evalJacobian().
50 */
51 void eval(double* x0, double* resid0, double rdt);
52
53 void reset() override;
54 const string type() const override { return "banded-direct"; }
55 void setValue(size_t row, size_t col, double value) override;
56 void initialize(size_t nVars) override;
57 void setBandwidth(size_t bw) override;
58 void updateTransient(double rdt, integer* mask) override;
59 void factorize() override {
60 m_mat.factor();
61 }
62
63 double& value(size_t i, size_t j) {
64 return m_mat.value(i, j);
65 }
66
67 double value(size_t i, size_t j) const {
68 return m_mat.value(i, j);
69 }
70
71 void solve(const double* const b, double* const x) {
72 m_mat.solve(b, x);
73 }
74
75 void solve(const size_t stateSize, double* b, double* x) override {
76 m_mat.solve(b, x);
77 }
78
79 int info() const override {
80 return m_mat.info();
81 }
82
83 //! Return the transient mask.
84 //! @deprecated Unused. To be removed after %Cantera 3.2.
85 vector<int>& transientMask() {
86 warn_deprecated("MultiJac::transientMask", "To be removed after Cantera 3.2");
87 return m_mask;
88 }
89
90protected:
91 //! Residual evaluator for this Jacobian
92 /*!
93 * This is a pointer to the residual evaluator. This object isn't owned by
94 * this Jacobian object.
95 *
96 * @deprecated Unused. To be removed after %Cantera 3.2.
97 */
98 OneDim* m_resid = nullptr;
99
100 BandMatrix m_mat; //!< Underlying matrix storage
101 vector<double> m_ssdiag; //!< Diagonal of the steady-state Jacobian
102
103 //! Transient mask for transient terms, 1 if transient, 0 if steady-state.
104 //! @deprecated Unused. To be removed after %Cantera 3.2.
105 vector<int> m_mask;
106};
107
108}
109
110#endif
Declarations for the class BandMatrix which is a child class of GeneralMatrix for banded matrices han...
Declarations for class SystemJacobian.
A class for banded matrices, involving matrix inversion processes.
Definition BandMatrix.h:37
int factor() override
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
int info() const
LAPACK "info" flag after last factor/solve operation.
Definition BandMatrix.h:239
int solve(const double *const b, double *const x)
Solve the matrix problem Ax = b.
double & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Definition MultiJac.h:25
BandMatrix m_mat
Underlying matrix storage.
Definition MultiJac.h:100
vector< int > & transientMask()
Return the transient mask.
Definition MultiJac.h:85
void setBandwidth(size_t bw) override
Used to provide system bandwidth for implementations that use banded matrix storage.
Definition MultiJac.cpp:35
void setValue(size_t row, size_t col, double value) override
Set a value at the specified row and column of the jacobian triplet vector.
Definition MultiJac.cpp:40
void factorize() override
Factorize the system matrix.
Definition MultiJac.h:59
vector< int > m_mask
Transient mask for transient terms, 1 if transient, 0 if steady-state.
Definition MultiJac.h:105
void eval(double *x0, double *resid0, double rdt)
Evaluates the Jacobian at x0 using finite differences.
Definition MultiJac.cpp:56
void solve(const size_t stateSize, double *b, double *x) override
Solve a linear system using the system matrix M
Definition MultiJac.h:75
vector< double > m_ssdiag
Diagonal of the steady-state Jacobian.
Definition MultiJac.h:101
const string type() const override
Derived type, corresponding to names registered with SystemJacobianFactory.
Definition MultiJac.h:54
OneDim * m_resid
Residual evaluator for this Jacobian.
Definition MultiJac.h:98
void reset() override
Reset parameters as needed.
Definition MultiJac.cpp:21
void updateTransient(double rdt, integer *mask) override
Update the diagonal terms in the Jacobian by using the transient mask .
Definition MultiJac.cpp:48
int info() const override
Get latest status of linear solver.
Definition MultiJac.h:79
void initialize(size_t nVars) override
Called during setup for any processes that need to be completed prior to setup functions used in sund...
Definition MultiJac.cpp:27
Container class for multiple-domain 1D problems.
Definition OneDim.h:28
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
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