svZeroDSolver
|
Sparse system. More...
#include <SparseSystem.h>
Public Member Functions | |
SparseSystem () | |
Construct a new Sparse System object. | |
SparseSystem (int n) | |
Construct a new Sparse System object. | |
~SparseSystem () | |
Destroy the Sparse System object. | |
void | reserve (Model *model) |
Reserve memory in system matrices based on number of triplets. | |
void | update_residual (Eigen::Matrix< double, Eigen::Dynamic, 1 > &y, Eigen::Matrix< double, Eigen::Dynamic, 1 > &ydot) |
Update the residual of the system. | |
void | update_jacobian (double time_coeff_ydot, double time_coeff_y) |
Update the jacobian of the system. | |
void | solve () |
Solve the system. | |
void | clean () |
Delete dynamically allocated memory (class member Eigen::SparseLU<Eigen::SparseMatrix> *solver) | |
Public Attributes | |
Eigen::SparseMatrix< double > | F |
System matrix F. | |
Eigen::SparseMatrix< double > | E |
System matrix E. | |
Eigen::SparseMatrix< double > | dC_dy |
System matrix dC/dy. | |
Eigen::SparseMatrix< double > | dC_dydot |
System matrix dC/dydot. | |
Eigen::Matrix< double, Eigen::Dynamic, 1 > | C |
System vector C. | |
Eigen::SparseMatrix< double > | jacobian |
Jacobian of the system. | |
Eigen::Matrix< double, Eigen::Dynamic, 1 > | residual |
Residual of the system. | |
Eigen::Matrix< double, Eigen::Dynamic, 1 > | dydot |
Solution increment of the system. | |
std::shared_ptr< Eigen::SparseLU< Eigen::SparseMatrix< double > > > | solver |
Sparse system.
This class contains all attributes and methods to create, modify, and solve sparse systems.
Flow rate, pressure, and other hemodynamic quantities in 0D models of vascular anatomies are governed by a system of nonlinear differential-algebraic equations (DAEs):
where and . Here, is the residual, is the vector of solution quantities and is its time derivative. is the total number of equations and the total number of global unknowns. The DAE system is solved implicitly using the generalized- method in Integrator. We then use the Newton-Raphson method to iteratively solve
with solution increment in iteration . The linearization of the time-discretized system is
with time factors and provided by Integrator.
SparseSystem::SparseSystem | ( | ) |
Construct a new Sparse System object.
SparseSystem::SparseSystem | ( | int | n | ) |
Construct a new Sparse System object.
n | Size of the system |
SparseSystem::~SparseSystem | ( | ) |
Destroy the Sparse System object.
void SparseSystem::clean | ( | ) |
Delete dynamically allocated memory (class member Eigen::SparseLU<Eigen::SparseMatrix> *solver)
void SparseSystem::reserve | ( | Model * | model | ) |
Reserve memory in system matrices based on number of triplets.
model | The model to reserve space for in the system |
void SparseSystem::solve | ( | ) |
Solve the system.
void SparseSystem::update_jacobian | ( | double | time_coeff_ydot, |
double | time_coeff_y ) |
Update the jacobian of the system.
time_coeff_ydot | Coefficent ydot-dependent part of jacobian |
time_coeff_y | Coefficent ydot-dependent part of jacobian |
void SparseSystem::update_residual | ( | Eigen::Matrix< double, Eigen::Dynamic, 1 > & | y, |
Eigen::Matrix< double, Eigen::Dynamic, 1 > & | ydot ) |
Update the residual of the system.
y | Vector of current solution quantities |
ydot | Derivate of y |
Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::C |
System vector C.
Eigen::SparseMatrix<double> SparseSystem::dC_dy |
System matrix dC/dy.
Eigen::SparseMatrix<double> SparseSystem::dC_dydot |
System matrix dC/dydot.
Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::dydot |
Solution increment of the system.
Eigen::SparseMatrix<double> SparseSystem::E |
System matrix E.
Eigen::SparseMatrix<double> SparseSystem::F |
System matrix F.
Eigen::SparseMatrix<double> SparseSystem::jacobian |
Jacobian of the system.
Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::residual |
Residual of the system.
std::shared_ptr<Eigen::SparseLU<Eigen::SparseMatrix<double> > > SparseSystem::solver |
Linear solver