Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplied by an instance of class OneDim. More...
#include <MultiJac.h>
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplied by an instance of class OneDim.
The residual function may consist of several linked 1D domains, with different variables in each domain.
Definition at line 23 of file MultiJac.h.
Public Member Functions | |
MultiJac (OneDim &r) | |
Constructor. | |
void | eval (double *x0, double *resid0, double rdt) |
Evaluates the Jacobian at x0 using finite differences. | |
double | elapsedTime () const |
Elapsed CPU time spent computing the Jacobian. | |
int | nEvals () const |
Number of Jacobian evaluations. | |
int | age () const |
Number of times 'incrementAge' has been called since the last evaluation. | |
void | incrementAge () |
Increment the Jacobian age. | |
void | updateTransient (double rdt, integer *mask) |
Update the transient terms in the Jacobian by using the transient mask. | |
void | setAge (int age) |
Set the Jacobian age. | |
vector< int > & | transientMask () |
Return the transient mask. | |
void | incrementDiagonal (int j, double d) |
Public Member Functions inherited from BandMatrix | |
BandMatrix () | |
Base Constructor. | |
BandMatrix (size_t n, size_t kl, size_t ku, double v=0.0) | |
Creates a banded matrix and sets all elements to zero. | |
BandMatrix (const BandMatrix &y) | |
BandMatrix & | operator= (const BandMatrix &y) |
void | resize (size_t n, size_t kl, size_t ku, double v=0.0) |
Resize the matrix problem. | |
void | bfill (double v=0.0) |
Fill or zero the matrix. | |
double & | operator() (size_t i, size_t j) override |
Index into the (i,j) element. | |
double | operator() (size_t i, size_t j) const override |
Constant Index into the (i,j) element. | |
double & | value (size_t i, size_t j) |
Return a changeable reference to element (i,j). | |
double | value (size_t i, size_t j) const |
Return the value of element (i,j). | |
size_t | index (size_t i, size_t j) const |
Returns the location in the internal 1D array corresponding to the (i,j) element in the banded array. | |
double | _value (size_t i, size_t j) const |
Return the value of the (i,j) element for (i,j) within the bandwidth. | |
size_t | nRows () const override |
Return the number of rows in the matrix. | |
size_t | nColumns () const |
Number of columns. | |
size_t | nSubDiagonals () const |
Number of subdiagonals. | |
size_t | nSuperDiagonals () const |
Number of superdiagonals. | |
size_t | ldim () const |
Return the number of rows of storage needed for the band storage. | |
void | mult (const double *b, double *prod) const override |
Multiply A*b and write result to prod . | |
void | leftMult (const double *const b, double *const prod) const override |
Multiply b*A and write result to prod. | |
int | factor () override |
Perform an LU decomposition, the LAPACK routine DGBTRF is used. | |
int | solve (const double *const b, double *const x) |
Solve the matrix problem Ax = b. | |
int | solve (double *b, size_t nrhs=1, size_t ldb=0) override |
Solve the matrix problem Ax = b. | |
void | zero () override |
Zero the matrix elements. | |
double | rcond (double a1norm) override |
Returns an estimate of the inverse of the condition number for the matrix. | |
int | factorAlgorithm () const override |
Returns the factor algorithm used. | |
double | oneNorm () const override |
Returns the one norm of the matrix. | |
double * | ptrColumn (size_t j) override |
Return a pointer to the top of column j. | |
double *const * | colPts () override |
Return a vector of const pointers to the columns. | |
size_t | checkRows (double &valueSmall) const override |
Check to see if we have any zero rows in the Jacobian. | |
size_t | checkColumns (double &valueSmall) const override |
Check to see if we have any zero columns in the Jacobian. | |
int | info () const |
LAPACK "info" flag after last factor/solve operation. | |
Public Member Functions inherited from GeneralMatrix | |
GeneralMatrix ()=default | |
Base Constructor. | |
virtual void | zero ()=0 |
Zero the matrix elements. | |
virtual void | mult (const double *b, double *prod) const =0 |
Multiply A*b and write result to prod. | |
virtual void | leftMult (const double *const b, double *const prod) const =0 |
Multiply b*A and write result to prod. | |
virtual int | factor ()=0 |
Factors the A matrix, overwriting A. | |
virtual int | factorQR () |
Factors the A matrix using the QR algorithm, overwriting A. | |
virtual double | rcondQR () |
Returns an estimate of the inverse of the condition number for the matrix. | |
virtual double | rcond (double a1norm)=0 |
Returns an estimate of the inverse of the condition number for the matrix. | |
virtual void | useFactorAlgorithm (int fAlgorithm) |
Change the way the matrix is factored. | |
virtual int | factorAlgorithm () const =0 |
Return the factor algorithm used. | |
virtual double | oneNorm () const =0 |
Calculate the one norm of the matrix. | |
virtual size_t | nRows () const =0 |
Return the number of rows in the matrix. | |
virtual void | clearFactorFlag () |
clear the factored flag | |
virtual int | solve (double *b, size_t nrhs=1, size_t ldb=0)=0 |
Solves the Ax = b system returning x in the b spot. | |
virtual bool | factored () const |
true if the current factorization is up to date with the matrix | |
virtual double * | ptrColumn (size_t j)=0 |
Return a pointer to the top of column j, columns are assumed to be contiguous in memory. | |
virtual double & | operator() (size_t i, size_t j)=0 |
Index into the (i,j) element. | |
virtual double | operator() (size_t i, size_t j) const =0 |
Constant Index into the (i,j) element. | |
virtual double *const * | colPts ()=0 |
Return a vector of const pointers to the columns. | |
virtual size_t | checkRows (double &valueSmall) const =0 |
Check to see if we have any zero rows in the Jacobian. | |
virtual size_t | checkColumns (double &valueSmall) const =0 |
Check to see if we have any zero columns in the Jacobian. | |
Protected Attributes | |
OneDim * | m_resid |
Residual evaluator for this Jacobian. | |
vector< double > | m_r1 |
Perturbed residual vector. | |
double | m_rtol = 1e-5 |
Relative tolerance for perturbing solution components. | |
double | m_atol = sqrt(std::numeric_limits<double>::epsilon()) |
Absolute tolerance for perturbing solution components. | |
double | m_elapsed = 0.0 |
Elapsed CPU time taken to compute the Jacobian. | |
vector< double > | m_ssdiag |
Diagonal of the steady-state Jacobian. | |
vector< int > | m_mask |
Transient mask for transient terms, 1 if transient, 0 if steady-state. | |
int | m_nevals = 0 |
Number of Jacobian evaluations. | |
int | m_age = 100000 |
Age of the Jacobian (times incrementAge() has been called) | |
Protected Attributes inherited from BandMatrix | |
vector< double > | data |
Matrix data. | |
vector< double > | ludata |
Factorized data. | |
size_t | m_n = 0 |
Number of rows and columns of the matrix. | |
size_t | m_kl = 0 |
Number of subdiagonals of the matrix. | |
size_t | m_ku = 0 |
Number of super diagonals of the matrix. | |
double | m_zero = 0 |
value of zero | |
unique_ptr< PivData > | m_ipiv |
Pivot vector. | |
vector< double * > | m_colPtrs |
Vector of column pointers. | |
vector< double * > | m_lu_col_ptrs |
vector< int > | iwork_ |
Extra work array needed - size = n. | |
vector< double > | work_ |
Extra dp work array needed - size = 3n. | |
int | m_info = 0 |
Protected Attributes inherited from GeneralMatrix | |
int | m_factored = false |
Indicates whether the matrix is factored. | |
Constructor.
r | The nonlinear system for which to compute the Jacobian. |
Definition at line 12 of file MultiJac.cpp.
void eval | ( | double * | x0, |
double * | resid0, | ||
double | rdt | ||
) |
Evaluates the Jacobian at x0 using finite differences.
The unperturbed residual function is resid0, which must be supplied on input. The parameter 'rdt' is the reciprocal of the time step. If zero, the steady-state Jacobian is evaluated, otherwise the transient Jacobian is evaluated.
Each variable in the input vector x0
is perturbed to compute the corresponding column of the Jacobian matrix. The Jacobian is computed by perturbing each variable, evaluating the residual function, and then estimating the partial derivatives numerically using finite differences.
x0 | Solution vector at which to evaluate the Jacobian |
resid0 | Residual vector at x0 |
rdt | Reciprocal of the time step |
Definition at line 35 of file MultiJac.cpp.
|
inline |
Elapsed CPU time spent computing the Jacobian.
Definition at line 48 of file MultiJac.h.
|
inline |
Number of Jacobian evaluations.
Definition at line 53 of file MultiJac.h.
|
inline |
Number of times 'incrementAge' has been called since the last evaluation.
Definition at line 58 of file MultiJac.h.
|
inline |
Increment the Jacobian age.
Definition at line 63 of file MultiJac.h.
void updateTransient | ( | double | rdt, |
integer * | mask | ||
) |
Update the transient terms in the Jacobian by using the transient mask.
Definition at line 21 of file MultiJac.cpp.
|
inline |
Set the Jacobian age.
Definition at line 71 of file MultiJac.h.
|
inline |
Return the transient mask.
Definition at line 76 of file MultiJac.h.
void incrementDiagonal | ( | int | j, |
double | d | ||
) |
Definition at line 28 of file MultiJac.cpp.
|
protected |
Residual evaluator for this Jacobian.
This is a pointer to the residual evaluator. This object isn't owned by this Jacobian object.
Definition at line 89 of file MultiJac.h.
|
protected |
Perturbed residual vector.
Definition at line 91 of file MultiJac.h.
|
protected |
Relative tolerance for perturbing solution components.
Definition at line 92 of file MultiJac.h.
|
protected |
Absolute tolerance for perturbing solution components.
Definition at line 95 of file MultiJac.h.
|
protected |
Elapsed CPU time taken to compute the Jacobian.
Definition at line 97 of file MultiJac.h.
|
protected |
Diagonal of the steady-state Jacobian.
Definition at line 98 of file MultiJac.h.
|
protected |
Transient mask for transient terms, 1 if transient, 0 if steady-state.
Definition at line 101 of file MultiJac.h.
|
protected |
Number of Jacobian evaluations.
Definition at line 103 of file MultiJac.h.
|
protected |
Age of the Jacobian (times incrementAge() has been called)
Definition at line 104 of file MultiJac.h.