System Jacobians that use Eigen sparse matrices for storage. More...
#include <EigenSparseJacobian.h>
System Jacobians that use Eigen sparse matrices for storage.
Definition at line 16 of file EigenSparseJacobian.h.
Public Member Functions | |
void | initialize (size_t networkSize) override |
Called during setup for any processes that need to be completed prior to setup functions used in sundials. | |
void | reset () override |
Reset parameters as needed. | |
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 | updatePreconditioner () override |
Transform Jacobian vector and write into preconditioner, P = (I - gamma * J) | |
void | updateTransient (double rdt, int *mask) override |
Update the diagonal terms in the Jacobian by using the transient mask \( \alpha \). | |
Eigen::SparseMatrix< double > | jacobian () |
Return underlying Jacobian matrix. | |
Eigen::SparseMatrix< double > | matrix () |
Return the internal preconditioner matrix. | |
void | printPreconditioner () override |
Print preconditioner contents. | |
void | printJacobian () |
Print jacobian contents. | |
![]() | |
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 Attributes | |
vector< Eigen::Triplet< double > > | m_jac_trips |
Vector of triples representing the jacobian used in preconditioning. | |
Eigen::SparseMatrix< double > | m_identity |
Storage of appropriately sized identity matrix for making the preconditioner. | |
Eigen::SparseMatrix< double > | m_matrix |
Container that is the sparse preconditioner. | |
![]() | |
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. | |
Additional Inherited Members | |
![]() | |
virtual void | factorize () |
Factorize the system matrix. | |
|
overridevirtual |
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 from SystemJacobian.
Definition at line 12 of file EigenSparseJacobian.cpp.
|
overridevirtual |
Reset parameters as needed.
Reimplemented from SystemJacobian.
Definition at line 30 of file EigenSparseJacobian.cpp.
|
overridevirtual |
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 from SystemJacobian.
Definition at line 36 of file EigenSparseJacobian.cpp.
|
overridevirtual |
Transform Jacobian vector and write into preconditioner, P = (I - gamma * J)
Reimplemented from SystemJacobian.
Definition at line 41 of file EigenSparseJacobian.cpp.
|
overridevirtual |
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 from SystemJacobian.
Definition at line 51 of file EigenSparseJacobian.cpp.
|
inline |
Return the internal preconditioner matrix.
Definition at line 31 of file EigenSparseJacobian.h.
|
overridevirtual |
Print preconditioner contents.
Reimplemented from SystemJacobian.
Definition at line 68 of file EigenSparseJacobian.cpp.
void printJacobian | ( | ) |
Print jacobian contents.
Definition at line 75 of file EigenSparseJacobian.cpp.
|
protected |
Vector of triples representing the jacobian used in preconditioning.
Definition at line 43 of file EigenSparseJacobian.h.
|
protected |
Storage of appropriately sized identity matrix for making the preconditioner.
Definition at line 46 of file EigenSparseJacobian.h.
|
protected |
Container that is the sparse preconditioner.
Definition at line 49 of file EigenSparseJacobian.h.