51 void eval(
double* x0,
double* resid0,
double rdt);
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;
63 double& value(
size_t i,
size_t j) {
67 double value(
size_t i,
size_t j)
const {
71 void solve(
const double*
const b,
double*
const x) {
75 void solve(
const size_t stateSize,
double* b,
double* x)
override {
79 int info()
const override {
86 warn_deprecated(
"MultiJac::transientMask",
"To be removed after Cantera 3.2");
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.
int factor() override
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
int info() const
LAPACK "info" flag after last factor/solve operation.
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...
BandMatrix m_mat
Underlying matrix storage.
vector< int > & transientMask()
Return the transient mask.
void setBandwidth(size_t bw) override
Used to provide system bandwidth for implementations that use banded matrix storage.
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.
void factorize() override
Factorize the system matrix.
vector< int > m_mask
Transient mask for transient terms, 1 if transient, 0 if steady-state.
void eval(double *x0, double *resid0, double rdt)
Evaluates the Jacobian at x0 using finite differences.
void solve(const size_t stateSize, double *b, double *x) override
Solve a linear system using the system matrix M
vector< double > m_ssdiag
Diagonal of the steady-state Jacobian.
const string type() const override
Derived type, corresponding to names registered with SystemJacobianFactory.
OneDim * m_resid
Residual evaluator for this Jacobian.
void reset() override
Reset parameters as needed.
void updateTransient(double rdt, integer *mask) override
Update the diagonal terms in the Jacobian by using the transient mask .
int info() const override
Get latest status of linear solver.
void initialize(size_t nVars) override
Called during setup for any processes that need to be completed prior to setup functions used in sund...
Container class for multiple-domain 1D problems.
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
Namespace for the Cantera kernel.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.