6#ifndef SYSTEMJACOBIAN_H
7#define SYSTEMJACOBIAN_H
51 virtual const string type()
const = 0;
57 virtual void setValue(
size_t row,
size_t col,
double value) {
102 virtual void solve(
const size_t stateSize,
double* rhs_vector,
double* output) {
233 "To be removed after Cantera 3.2. Renamed to SystemJacobian.");
An error indicating that an unimplemented function has been called.
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
virtual double gamma()
Get gamma used in preconditioning.
int nEvals() const
Number of Jacobian evaluations.
void setAge(int age)
Set the Jacobian age.
virtual void setGamma(double gamma)
Set gamma used in preconditioning.
size_t m_dim
Dimension of the system.
void incrementEvals()
Increment the number of times the Jacobian has been evaluated.
string preconditionerSide() const
Get preconditioner application side for CVODES.
virtual void setup()
Perform preconditioner specific post-reactor setup operations such as factorize.
double m_elapsed
Elapsed CPU time taken to compute the Jacobian.
double m_gamma
gamma value used in M = I - gamma*J
virtual int info() const
Get latest status of linear solver.
virtual void setPreconditionerSide(const string &preconSide)
For iterative solvers, set the side where the preconditioner is applied.
virtual void printPreconditioner()
Print preconditioner contents.
virtual void factorize()
Factorize the system matrix.
virtual void solve(const size_t stateSize, double *rhs_vector, double *output)
Solve a linear system using the system matrix M
virtual void reset()
Reset parameters as needed.
string m_precon_side
For iterative solvers, side of the system to apply the preconditioner on.
virtual void setBandwidth(size_t bw)
Used to provide system bandwidth for implementations that use banded matrix storage.
void incrementAge()
Increment the Jacobian age.
virtual void stateAdjustment(vector< double > &state)
Adjust the state vector based on the preconditioner, e.g., Adaptive preconditioning uses a strictly p...
void updateElapsed(double evalTime)
Increase the elapsed time by the specified time.
int age() const
Number of times 'incrementAge' has been called since the last evaluation.
virtual const string type() const =0
Derived type, corresponding to names registered with SystemJacobianFactory.
int m_nevals
Number of Jacobian evaluations.
bool m_init
bool saying whether or not the system is initialized
int m_age
Age of the Jacobian (times incrementAge() has been called)
void clearStats()
Clear collected stats about number of Jacobian evaluations, CPU time spent on Jacobian updates,...
virtual void initialize(size_t networkSize)
Called during setup for any processes that need to be completed prior to setup functions used in sund...
double m_atol
Absolute tolerance of the ODE solver.
virtual void setAbsoluteTolerance(double atol)
Set the absolute tolerance in the solver outside of the network initialization.
virtual void setValue(size_t row, size_t col, double value)
Set a value at the specified row and column of the jacobian triplet vector.
double elapsedTime() const
Elapsed CPU time spent computing the Jacobian elements.
virtual void updatePreconditioner()
Transform Jacobian vector and write into preconditioner, P = (I - gamma * J)
virtual void updateTransient(double rdt, int *mask)
Update the diagonal terms in the Jacobian by using the transient mask .
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
Namespace for the Cantera kernel.
PreconditionerSide
Specifies the side of the system on which the preconditioner is applied.
@ LEFT_PRECONDITION
No preconditioning.
@ BOTH_PRECONDITION
Right side preconditioning.
@ RIGHT_PRECONDITION
Left side preconditioning.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.