Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers. More...
#include <SystemJacobian.h>
Abstract base class representing Jacobian matrices and preconditioners used in nonlinear solvers.
The matrix is initially populated with elements of the (steady-state) Jacobian using the setValue() method.
If the object is being used for preconditioning, the preconditioner matrix M = I - \gamma J is computed and factorized by calling updatePreconditioner(). If the object is being used for solving a (pseudo-)transient problem, the transient Jacobian M = J - \alpha \Delta t^{-1} is computed and factorized by calling the updateTransient() method.
Once the system is factorized, the solve() method can be used for preconditioning or solving Newton steps as appropriate.
Implementations of this class provide different options for how to store the matrix elements and different algorithms for factorizing the matrix.
Definition at line 43 of file SystemJacobian.h.
Public Member Functions | |
virtual const string | type () const =0 |
Derived type, corresponding to names registered with SystemJacobianFactory. | |
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. | |
virtual void | stateAdjustment (vector< double > &state) |
Adjust the state vector based on the preconditioner, e.g., Adaptive preconditioning uses a strictly positive composition when preconditioning which is handled by this function. | |
string | preconditionerSide () const |
Get preconditioner application side for CVODES. | |
virtual void | setPreconditionerSide (const string &preconSide) |
For iterative solvers, set the side where the preconditioner is applied. | |
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 \alpha . | |
virtual void | solve (const size_t stateSize, double *rhs_vector, double *output) |
Solve a linear system using the system matrix M | |
virtual void | setup () |
Perform preconditioner specific post-reactor setup operations such as factorize. | |
virtual void | reset () |
Reset parameters as needed. | |
virtual void | initialize (size_t networkSize) |
Called during setup for any processes that need to be completed prior to setup functions used in sundials. | |
virtual void | setBandwidth (size_t bw) |
Used to provide system bandwidth for implementations that use banded matrix storage. | |
virtual void | printPreconditioner () |
Print preconditioner contents. | |
virtual void | setGamma (double gamma) |
Set gamma used in preconditioning. | |
virtual double | gamma () |
Get gamma used in preconditioning. | |
virtual void | setAbsoluteTolerance (double atol) |
Set the absolute tolerance in the solver outside of the network initialization. | |
virtual int | info () const |
Get latest status of linear solver. | |
double | elapsedTime () const |
Elapsed CPU time spent computing the Jacobian elements. | |
void | updateElapsed (double evalTime) |
Increase the elapsed time by the specified time. | |
int | nEvals () const |
Number of Jacobian evaluations. | |
void | incrementEvals () |
Increment the number of times the Jacobian has been evaluated. | |
int | age () const |
Number of times 'incrementAge' has been called since the last evaluation. | |
void | incrementAge () |
Increment the Jacobian age. | |
void | setAge (int age) |
Set the Jacobian age. | |
void | clearStats () |
Clear collected stats about number of Jacobian evaluations, CPU time spent on Jacobian updates, and the number of Newton steps that have been taken since updating the Jacobian. | |
Protected Member Functions | |
virtual void | factorize () |
Factorize the system matrix. | |
Protected Attributes | |
size_t | m_dim |
Dimension of the system. | |
double | m_gamma = 1.0 |
gamma value used in M = I - gamma*J | |
bool | m_init = false |
bool saying whether or not the system is initialized | |
double | m_atol = 0 |
Absolute tolerance of the ODE solver. | |
double | m_elapsed = 0.0 |
Elapsed CPU time taken to compute the Jacobian. | |
int | m_nevals = 0 |
Number of Jacobian evaluations. | |
int | m_age = 100000 |
Age of the Jacobian (times incrementAge() has been called) | |
string | m_precon_side = "none" |
For iterative solvers, side of the system to apply the preconditioner on. | |
|
inline |
Definition at line 46 of file SystemJacobian.h.
|
inlinevirtual |
Definition at line 48 of file SystemJacobian.h.
|
pure virtual |
Derived type, corresponding to names registered with SystemJacobianFactory.
Implemented in AdaptivePreconditioner, EigenSparseDirectJacobian, and MultiJac.
|
inlinevirtual |
Set a value at the specified row and column of the jacobian triplet vector.
row | row in the jacobian matrix |
col | column in the jacobian matrix |
value | value of the element at row and col |
Reimplemented in EigenSparseJacobian, and MultiJac.
Definition at line 57 of file SystemJacobian.h.
|
inlinevirtual |
Adjust the state vector based on the preconditioner, e.g., Adaptive preconditioning uses a strictly positive composition when preconditioning which is handled by this function.
state | a vector containing the state to be updated |
Reimplemented in AdaptivePreconditioner.
Definition at line 65 of file SystemJacobian.h.
|
inline |
Get preconditioner application side for CVODES.
Definition at line 70 of file SystemJacobian.h.
|
inlinevirtual |
For iterative solvers, set the side where the preconditioner is applied.
Definition at line 75 of file SystemJacobian.h.
|
inlinevirtual |
Transform Jacobian vector and write into preconditioner, P = (I - gamma * J)
Reimplemented in EigenSparseJacobian.
Definition at line 81 of file SystemJacobian.h.
|
inlinevirtual |
Update the diagonal terms in the Jacobian by using the transient mask \alpha .
rdt | Reciprocal of the time step [1/s] |
mask | Mask for transient terms: 1 if transient, 0 if algebraic. |
Reimplemented in EigenSparseJacobian, and MultiJac.
Definition at line 89 of file SystemJacobian.h.
|
inlinevirtual |
Solve a linear system using the system matrix M
The matrix M can either be the preconditioner or the transient Jacobian matrix, depending on whether updateTransient() or updatePreconditioner() was called last.
[in] | stateSize | length of the rhs and output vectors |
[in] | rhs_vector | right hand side vector used in linear system |
[out] | output | output vector for solution |
Reimplemented in MultiJac, AdaptivePreconditioner, and EigenSparseDirectJacobian.
Definition at line 102 of file SystemJacobian.h.
|
inlinevirtual |
Perform preconditioner specific post-reactor setup operations such as factorize.
Reimplemented in AdaptivePreconditioner.
Definition at line 109 of file SystemJacobian.h.
|
inlinevirtual |
Reset parameters as needed.
Reimplemented in EigenSparseJacobian, and MultiJac.
Definition at line 114 of file SystemJacobian.h.
|
inlinevirtual |
Called during setup for any processes that need to be completed prior to setup functions used in sundials.
networkSize | the number of variables in the associated reactor network |
Reimplemented in AdaptivePreconditioner, EigenSparseJacobian, and MultiJac.
Definition at line 121 of file SystemJacobian.h.
|
inlinevirtual |
Used to provide system bandwidth for implementations that use banded matrix storage.
Ignored if not needed.
Reimplemented in MultiJac.
Definition at line 127 of file SystemJacobian.h.
|
inlinevirtual |
Print preconditioner contents.
Reimplemented in EigenSparseJacobian.
Definition at line 130 of file SystemJacobian.h.
|
inlinevirtual |
Set gamma used in preconditioning.
gamma | used in M = I - gamma*J |
Definition at line 136 of file SystemJacobian.h.
|
inlinevirtual |
Get gamma used in preconditioning.
Definition at line 141 of file SystemJacobian.h.
|
inlinevirtual |
Set the absolute tolerance in the solver outside of the network initialization.
atol | the specified tolerance |
Definition at line 147 of file SystemJacobian.h.
|
inlinevirtual |
Get latest status of linear solver.
Zero indicates success. Meaning of non-zero values is solver dependent.
Reimplemented in AdaptivePreconditioner, and MultiJac.
Definition at line 153 of file SystemJacobian.h.
|
inline |
Elapsed CPU time spent computing the Jacobian elements.
Definition at line 158 of file SystemJacobian.h.
|
inline |
Increase the elapsed time by the specified time.
Definition at line 163 of file SystemJacobian.h.
|
inline |
Number of Jacobian evaluations.
Definition at line 168 of file SystemJacobian.h.
|
inline |
Increment the number of times the Jacobian has been evaluated.
Definition at line 172 of file SystemJacobian.h.
|
inline |
Number of times 'incrementAge' has been called since the last evaluation.
Definition at line 177 of file SystemJacobian.h.
|
inline |
Increment the Jacobian age.
Definition at line 182 of file SystemJacobian.h.
|
inline |
Set the Jacobian age.
Definition at line 187 of file SystemJacobian.h.
|
inline |
Clear collected stats about number of Jacobian evaluations, CPU time spent on Jacobian updates, and the number of Newton steps that have been taken since updating the Jacobian.
Definition at line 194 of file SystemJacobian.h.
|
inlineprotectedvirtual |
Factorize the system matrix.
This method should be called at the end of implementations of updatePreconditioner() and updateTransient().
Reimplemented in AdaptivePreconditioner, EigenSparseDirectJacobian, and MultiJac.
Definition at line 203 of file SystemJacobian.h.
|
protected |
Dimension of the system.
Definition at line 208 of file SystemJacobian.h.
|
protected |
gamma value used in M = I - gamma*J
Definition at line 211 of file SystemJacobian.h.
|
protected |
bool saying whether or not the system is initialized
Definition at line 214 of file SystemJacobian.h.
|
protected |
Absolute tolerance of the ODE solver.
Definition at line 217 of file SystemJacobian.h.
|
protected |
Elapsed CPU time taken to compute the Jacobian.
Definition at line 219 of file SystemJacobian.h.
|
protected |
Number of Jacobian evaluations.
Definition at line 220 of file SystemJacobian.h.
|
protected |
Age of the Jacobian (times incrementAge() has been called)
Definition at line 222 of file SystemJacobian.h.
|
protected |
For iterative solvers, side of the system to apply the preconditioner on.
Definition at line 225 of file SystemJacobian.h.