svZeroDSolver
Loading...
Searching...
No Matches
SparseSystem Class Reference

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
 

Detailed Description

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):

\[
\mathbf{r}(\boldsymbol{\alpha}, \mathbf{y},\dot{\mathbf{y}}, t) =
\mathbf{E}(\boldsymbol{\alpha}) \cdot \dot{\mathbf{y}} +
\mathbf{F}(\boldsymbol{\alpha}) \cdot \mathbf{y} +
\mathbf{c}(\mathbf{y},\dot{\mathbf{y}}, t) = \mathbf{0}
\]

where $\mathbf{r},\textbf{y},\textbf{c} \in \mathbb{R}^{N}$ and $\textbf{E},\textbf{F} \in \mathbb{R}^{N \times N}$. Here, $\textbf{r}$ is the residual, $\textbf{y}$ is the vector of solution quantities and $\dot{\textbf{y}}$ is its time derivative. $N$ is the total number of equations and the total number of global unknowns. The DAE system is solved implicitly using the generalized- $\alpha$ method in Integrator. We then use the Newton-Raphson method to iteratively solve

\[
\mathbf{K}^{i} \cdot \Delta\dot{\mathbf{y}}^{i} = - \mathbf{r}^{i},
\]

with solution increment $\Delta\dot{\mathbf{y}}^{i}$ in iteration $i$. The linearization of the time-discretized system is

\[
\mathbf{K} =
\frac{\partial \mathbf{r}}{\partial \mathbf{y}} =
c_{\dot{\mathbf{y}}} \left( \mathbf{E} + \frac{\partial \mathbf{c}}{\partial
\dot{\mathbf{y}}} \right) +
c_{\mathbf{y}} \left( \mathbf{F} + \frac{\partial \mathbf{c}}{\partial
\mathbf{y}} \right), \]

with time factors $c_{\dot{\mathbf{y}}}=\alpha_m$ and $c_{\mathbf{y}}=\alpha_f\gamma\Delta t$ provided by Integrator.

Constructor & Destructor Documentation

◆ SparseSystem() [1/2]

SparseSystem::SparseSystem ( )

Construct a new Sparse System object.

◆ SparseSystem() [2/2]

SparseSystem::SparseSystem ( int n)

Construct a new Sparse System object.

Parameters
nSize of the system

◆ ~SparseSystem()

SparseSystem::~SparseSystem ( )

Destroy the Sparse System object.

Member Function Documentation

◆ clean()

void SparseSystem::clean ( )

Delete dynamically allocated memory (class member Eigen::SparseLU<Eigen::SparseMatrix> *solver)

◆ reserve()

void SparseSystem::reserve ( Model * model)

Reserve memory in system matrices based on number of triplets.

Parameters
modelThe model to reserve space for in the system

◆ solve()

void SparseSystem::solve ( )

Solve the system.

◆ update_jacobian()

void SparseSystem::update_jacobian ( double time_coeff_ydot,
double time_coeff_y )

Update the jacobian of the system.

Parameters
time_coeff_ydotCoefficent ydot-dependent part of jacobian
time_coeff_yCoefficent ydot-dependent part of jacobian

◆ update_residual()

void SparseSystem::update_residual ( Eigen::Matrix< double, Eigen::Dynamic, 1 > & y,
Eigen::Matrix< double, Eigen::Dynamic, 1 > & ydot )

Update the residual of the system.

Parameters
yVector of current solution quantities
ydotDerivate of y

Member Data Documentation

◆ C

Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::C

System vector C.

◆ dC_dy

Eigen::SparseMatrix<double> SparseSystem::dC_dy

System matrix dC/dy.

◆ dC_dydot

Eigen::SparseMatrix<double> SparseSystem::dC_dydot

System matrix dC/dydot.

◆ dydot

Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::dydot

Solution increment of the system.

◆ E

Eigen::SparseMatrix<double> SparseSystem::E

System matrix E.

◆ F

Eigen::SparseMatrix<double> SparseSystem::F

System matrix F.

◆ jacobian

Eigen::SparseMatrix<double> SparseSystem::jacobian

Jacobian of the system.

◆ residual

Eigen::Matrix<double, Eigen::Dynamic, 1> SparseSystem::residual

Residual of the system.

◆ solver

std::shared_ptr<Eigen::SparseLU<Eigen::SparseMatrix<double> > > SparseSystem::solver
Initial value:
=
std::shared_ptr<Eigen::SparseLU<Eigen::SparseMatrix<double>>>(
new Eigen::SparseLU<Eigen::SparseMatrix<double>>())

Linear solver


The documentation for this class was generated from the following files: