Cantera  3.3.0a1
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
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 MultiJac() = default;
28 void reset() override;
29 const string type() const override { return "banded-direct"; }
30 void setValue(size_t row, size_t col, double value) override;
31 void initialize(size_t nVars) override;
32 void setBandwidth(size_t bw) override;
33 void updateTransient(double rdt, integer* mask) override;
34 void factorize() override {
35 m_mat.factor();
36 }
37
38 double& value(size_t i, size_t j) {
39 return m_mat.value(i, j);
40 }
41
42 double value(size_t i, size_t j) const {
43 return m_mat.value(i, j);
44 }
45
46 void solve(const double* const b, double* const x) {
47 m_mat.solve(b, x);
48 }
49
50 void solve(const size_t stateSize, double* b, double* x) override {
51 m_mat.solve(b, x);
52 }
53
54 int info() const override {
55 return m_mat.info();
56 }
57
58protected:
59 BandMatrix m_mat; //!< Underlying matrix storage
60 vector<double> m_ssdiag; //!< Diagonal of the steady-state Jacobian
61};
62
63}
64
65#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
void solve(const double *const b, double *const x)
Solve the matrix problem Ax = b.
int info() const
LAPACK "info" flag after last factor/solve operation.
Definition BandMatrix.h:230
void factor() override
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
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:59
void setBandwidth(size_t bw) override
Used to provide system bandwidth for implementations that use banded matrix storage.
Definition MultiJac.cpp:25
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:30
void factorize() override
Factorize the system matrix.
Definition MultiJac.h:34
void solve(const size_t stateSize, double *b, double *x) override
Solve a linear system using the system matrix M
Definition MultiJac.h:50
vector< double > m_ssdiag
Diagonal of the steady-state Jacobian.
Definition MultiJac.h:60
const string type() const override
Derived type, corresponding to names registered with SystemJacobianFactory.
Definition MultiJac.h:29
void reset() override
Reset parameters as needed.
Definition MultiJac.cpp:12
void updateTransient(double rdt, integer *mask) override
Update the diagonal terms in the Jacobian by using the transient mask .
Definition MultiJac.cpp:38
int info() const override
Get latest status of linear solver.
Definition MultiJac.h:54
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:18
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595