|
svZeroDSolver
|
Closed loop coronary boundary condition which is connected to other blocks on both sides and the intramyocardial pressure is specified by the pressure in a heart block (not as a parameter). More...
#include <ClosedLoopCoronaryBC.h>
Public Types | |
| enum | ParamId { RA = 0 , RAM = 1 , RV = 2 , CA = 3 , CIM = 4 } |
| Local IDs of the parameters. More... | |
Public Member Functions | |
| ClosedLoopCoronaryBC (int id, Model *model, BlockType block_type) | |
| Construct a ClosedLoopCoronaryBC object. | |
| void | setup_dofs (DOFHandler &dofhandler) |
| Set up the degrees of freedom (DOF) of the block. | |
| virtual void | setup_model_dependent_params ()=0 |
| Setup parameters that depend on the model. | |
| void | update_constant (SparseSystem &system, std::vector< double > ¶meters) |
| Update the constant contributions of the element in a sparse system. | |
| void | update_solution (SparseSystem &system, std::vector< double > ¶meters, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &y, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &dy) |
| Update the solution-dependent contributions of the element in a sparse system. | |
Public Member Functions inherited from Block | |
| Block (int id, Model *model, BlockType block_type, BlockClass block_class, std::vector< std::pair< std::string, InputParameter > > input_params) | |
| Construct a new Block object. | |
| ~Block () | |
| Destroy the Block object. | |
| Block (const Block &)=delete | |
| Copy the Block object. | |
| std::string | get_name () |
| Get the name of the block. | |
| void | update_vessel_type (VesselType type) |
| Update vessel type of the block. | |
| void | setup_params_ (const std::vector< int > ¶m_ids) |
| Setup parameter IDs for the block. | |
| void | setup_dofs_ (DOFHandler &dofhandler, int num_equations, const std::list< std::string > &internal_var_names) |
| Set up the degrees of freedom (DOF) of the block. | |
| virtual void | setup_initial_state_dependent_params (State initial_state, std::vector< double > ¶meters) |
| Setup parameters that depend on the initial state. | |
| virtual void | update_time (SparseSystem &system, std::vector< double > ¶meters) |
| Update the time-dependent contributions of the element in a sparse system. | |
| virtual void | post_solve (Eigen::Matrix< double, Eigen::Dynamic, 1 > &y) |
| Modify the solution after solving it. | |
| virtual void | update_gradient (Eigen::SparseMatrix< double > &jacobian, Eigen::Matrix< double, Eigen::Dynamic, 1 > &residual, Eigen::Matrix< double, Eigen::Dynamic, 1 > &alpha, std::vector< double > &y, std::vector< double > &dy) |
| Set the gradient of the block contributions with respect to the parameters. | |
| virtual TripletsContributions | get_num_triplets () |
| Get number of triplets of element. | |
Public Attributes | |
| TripletsContributions | num_triplets {9, 5, 0} |
| Number of triplets of element. | |
Public Attributes inherited from Block | |
| const int | id |
| Global ID of the block. | |
| const Model * | model |
| The model to which the block belongs. | |
| const BlockType | block_type |
| Type of this block. | |
| const BlockClass | block_class |
| Class of this block. | |
| VesselType | vessel_type = VesselType::neither |
| Vessel type of this block. | |
| const std::vector< std::pair< std::string, InputParameter > > | input_params |
| Map from name to input parameter. | |
| std::vector< Node * > | inlet_nodes |
| Inlet nodes. | |
| std::vector< Node * > | outlet_nodes |
| Outlet nodes. | |
| bool | steady = false |
| Toggle steady behavior. | |
| bool | input_params_list = false |
| Are input parameters given as a list? | |
| std::vector< int > | global_param_ids |
| Global IDs for the block parameters. | |
| std::vector< int > | global_var_ids |
| Global variable indices of the local element contributions. | |
| std::vector< int > | global_eqn_ids |
| Global equation indices of the local element contributions. | |
| TripletsContributions | num_triplets |
| Number of triplets of element. | |
Protected Attributes | |
| int | ventricle_var_id |
| Variable index of either left or right ventricle. | |
| int | im_param_id |
| Index of parameter Im. | |
Closed loop coronary boundary condition which is connected to other blocks on both sides and the intramyocardial pressure is specified by the pressure in a heart block (not as a parameter).
![\[\begin{circuitikz} \draw
node[left] {$Q_{in}$} [-latex] (0,0) -- (0.8,0);
\draw (1,0) node[anchor=south]{$P_{in}$}
to [R, l=$R_a$, *-] (3,0)
to [R, l=$R_{am}$, -] (5,0)
to [R, l=$R_v$, *-*] (7,0)
node[anchor=south]{$P_{out}$}
(5,0) to [C, l=$C_{im} \;V_{im}$, -*] (5,-1.5)
node[left]{$P_{im}$}
(3,0) to [C, l=$C_a$, -*] (3,-1.5)
node[left]{$P_a$};
\draw [-latex] (7.2,0) -- (8.0,0) node[right] {$Q_{out}$};
\end{circuitikz}
\]](form_68.png)
![\[P_{out} - P_{in} + (R_{am}+R_a)Q_{in} + R_v Q_{out} + R_{am} C_a
\frac{dP_a}{dt} - R_{am} C_a \frac{dP_{in}}{dt} + R_{am} R_a C_a
\frac{dQ_{in}}{dt} = 0 \]](form_69.png)
![\[Q_{in} - Q_{out} + C_a \frac{dP_a}{dt} - C_a \frac{dP_{in}}{dt} + C_a R_a
\frac{dQ_{in}}{dt} - \frac{dV_{im}}{dt} = 0 \]](form_70.png)
![\[C_{im} P_{out} + C_{im} R_v Q_{out} - C_{im} P_{im} - V_{im} = 0
\]](form_71.png)
![\[\mathbf{y}^{e}=\left[\begin{array}{lllll}P_{in} & Q_{in} & P_{out} &
Q_{out} & V_{im}\end{array}\right]^{T}, \]](form_72.png)
![\[\mathbf{E}^{e}=\left[\begin{array}{ccccc}
-R_{am} C_{a} & R_{am} R_{a} C_{a} & 0 & 0 & 0 \\
-C_{a} & R_{a} C_{a} & 0 & 0 & -1 \\
0 & 0 & 0 & 0 & 0 \\
\end{array}\right] \]](form_73.png)
![\[\mathbf{F}^{e}=\left[\begin{array}{ccccc}
-1 & R_{am} + R_{a} & 1 & R_v & 0 \\
0 & 1 & 0 & -1 & 0 \\
0 & 0 & C_{im} & C_{im} R_v & -1 \\
\end{array}\right] \]](form_74.png)
![\[\mathbf{c}^{e}=\left[\begin{array}{c}
C_{a} R_{am} \frac{d P_{a}}{d t} \\
C_{a}\frac{d P_{a}}{d t} \\
-C_{im} P_{im}
\end{array}\right] \]](form_75.png)
Assume 
Parameter sequence for constructing this block
0 Ra: Small artery resistance1 Ram: Microvascular resistance2 Rv: Venous resistance3 Ca: Small artery capacitance4 Cim: Intramyocardial capacitanceNames of internal variables in this block's output:
volume_im: Intramyocardial volume Local IDs of the parameters.
Construct a ClosedLoopCoronaryBC object.
| id | Global ID of the block |
| model | The model to which the block belongs |
| block_type | The specific type of block (left or right) |
|
virtual |
Set up the degrees of freedom (DOF) of the block.
Set global_var_ids and global_eqn_ids of the element based on the number of equations and the number of internal variables of the element.
| dofhandler | Degree-of-freedom handler to register variables and equations at |
Reimplemented from Block.
|
pure virtual |
Setup parameters that depend on the model.
Reimplemented from Block.
Implemented in ClosedLoopCoronaryLeftBC, and ClosedLoopCoronaryRightBC.
|
virtual |
Update the constant contributions of the element in a sparse system.
| system | System to update contributions at |
| parameters | Parameters of the model |
Reimplemented from Block.
|
virtual |
Update the solution-dependent contributions of the element in a sparse system.
| system | System to update contributions at |
| parameters | Parameters of the model |
| y | Current solution |
| dy | Current derivate of the solution |
Reimplemented from Block.
|
protected |
Index of parameter Im.
| TripletsContributions ClosedLoopCoronaryBC::num_triplets {9, 5, 0} |
Number of triplets of element.
Number of triplets that the element contributes to the global system (relevant for sparse memory reservation)
|
protected |
Variable index of either left or right ventricle.