svZeroDSolver
Loading...
Searching...
No Matches
ClosedLoopCoronaryBC.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Stanford University, The Regents of the
2// University of California, and others. SPDX-License-Identifier: BSD-3-Clause
3/**
4 * @file ClosedLoopCoronaryBC.h
5 * @brief model::ClosedLoopCoronaryBC source file
6 */
7#ifndef SVZERODSOLVER_MODEL_CLOSEDLOOPCORONARYBC_HPP_
8#define SVZERODSOLVER_MODEL_CLOSEDLOOPCORONARYBC_HPP_
9
10#include "Block.h"
12#include "SparseSystem.h"
13
14/**
15 * @brief Closed loop coronary boundary condition which is connected to
16 * other blocks on both sides and the intramyocardial pressure is
17 * specified by the pressure in a heart block (not as a parameter).
18 *
19 * \f[
20 * \begin{circuitikz} \draw
21 * node[left] {$Q_{in}$} [-latex] (0,0) -- (0.8,0);
22 * \draw (1,0) node[anchor=south]{$P_{in}$}
23 * to [R, l=$R_a$, *-] (3,0)
24 * to [R, l=$R_{am}$, -] (5,0)
25 * to [R, l=$R_v$, *-*] (7,0)
26 * node[anchor=south]{$P_{out}$}
27 * (5,0) to [C, l=$C_{im} \;V_{im}$, -*] (5,-1.5)
28 * node[left]{$P_{im}$}
29 * (3,0) to [C, l=$C_a$, -*] (3,-1.5)
30 * node[left]{$P_a$};
31 * \draw [-latex] (7.2,0) -- (8.0,0) node[right] {$Q_{out}$};
32 * \end{circuitikz}
33 * \f]
34 *
35 * ### Governing equations
36 *
37 * \f[
38 * P_{out} - P_{in} + (R_{am}+R_a)Q_{in} + R_v Q_{out} + R_{am} C_a
39 * \frac{dP_a}{dt} - R_{am} C_a \frac{dP_{in}}{dt} + R_{am} R_a C_a
40 * \frac{dQ_{in}}{dt} = 0 \f]
41 *
42 * \f[
43 * Q_{in} - Q_{out} + C_a \frac{dP_a}{dt} - C_a \frac{dP_{in}}{dt} + C_a R_a
44 * \frac{dQ_{in}}{dt} - \frac{dV_{im}}{dt} = 0 \f]
45 *
46 * \f[
47 * C_{im} P_{out} + C_{im} R_v Q_{out} - C_{im} P_{im} - V_{im} = 0
48 * \f]
49 *
50 * ### Local contributions
51 *
52 * \f[
53 * \mathbf{y}^{e}=\left[\begin{array}{lllll}P_{in} & Q_{in} & P_{out} &
54 * Q_{out} & V_{im}\end{array}\right]^{T}, \f]
55 *
56 * \f[
57 * \mathbf{E}^{e}=\left[\begin{array}{ccccc}
58 * -R_{am} C_{a} & R_{am} R_{a} C_{a} & 0 & 0 & 0 \\
59 * -C_{a} & R_{a} C_{a} & 0 & 0 & -1 \\
60 * 0 & 0 & 0 & 0 & 0 \\
61 * \end{array}\right] \f]
62 *
63 *
64 * \f[
65 * \mathbf{F}^{e}=\left[\begin{array}{ccccc}
66 * -1 & R_{am} + R_{a} & 1 & R_v & 0 \\
67 * 0 & 1 & 0 & -1 & 0 \\
68 * 0 & 0 & C_{im} & C_{im} R_v & -1 \\
69 * \end{array}\right] \f]
70 *
71 * \f[
72 * \mathbf{c}^{e}=\left[\begin{array}{c}
73 * C_{a} R_{am} \frac{d P_{a}}{d t} \\
74 * C_{a}\frac{d P_{a}}{d t} \\
75 * -C_{im} P_{im}
76 * \end{array}\right] \f]
77 *
78 * Assume \f$P_a=0\f$.
79 *
80 * ### Parameters
81 *
82 * Parameter sequence for constructing this block
83 *
84 * * `0` Ra: Small artery resistance
85 * * `1` Ram: Microvascular resistance
86 * * `2` Rv: Venous resistance
87 * * `3` Ca: Small artery capacitance
88 * * `4` Cim: Intramyocardial capacitance
89 *
90 * ### Internal variables
91 *
92 * Names of internal variables in this block's output:
93 *
94 * * `volume_im`: Intramyocardial volume
95 *
96 */
98 public:
99 /**
100 * @brief Construct a ClosedLoopCoronaryBC object.
101 *
102 * @param id Global ID of the block
103 * @param model The model to which the block belongs
104 * @param block_type The specific type of block (left or right)
105 */
107 : Block(id, model, block_type, BlockClass::closed_loop,
108 {{"Ra", InputParameter()},
109 {"Ram", InputParameter()},
110 {"Rv", InputParameter()},
111 {"Ca", InputParameter()},
112 {"Cim", InputParameter()}}) {}
113
114 /**
115 * @brief Local IDs of the parameters
116 *
117 */
118 enum ParamId { RA = 0, RAM = 1, RV = 2, CA = 3, CIM = 4 };
119
120 /**
121 * @brief Set up the degrees of freedom (DOF) of the block
122 *
123 * Set \ref global_var_ids and \ref global_eqn_ids of the element based on the
124 * number of equations and the number of internal variables of the
125 * element.
126 *
127 * @param dofhandler Degree-of-freedom handler to register variables and
128 * equations at
129 */
130 void setup_dofs(DOFHandler &dofhandler);
131
132 /**
133 * @brief Setup parameters that depend on the model
134 *
135 */
137
138 /**
139 * @brief Update the constant contributions of the element in a sparse system
140 *
141 * @param system System to update contributions at
142 * @param parameters Parameters of the model
143 */
144 void update_constant(SparseSystem &system, std::vector<double> &parameters);
145
146 /**
147 * @brief Update the solution-dependent contributions of the element in a
148 * sparse system
149 *
150 * @param system System to update contributions at
151 * @param parameters Parameters of the model
152 * @param y Current solution
153 * @param dy Current derivate of the solution
154 */
155 void update_solution(SparseSystem &system, std::vector<double> &parameters,
156 const Eigen::Matrix<double, Eigen::Dynamic, 1> &y,
157 const Eigen::Matrix<double, Eigen::Dynamic, 1> &dy);
158
159 /**
160 * @brief Number of triplets of element
161 *
162 * Number of triplets that the element contributes to the global system
163 * (relevant for sparse memory reservation)
164 */
166
167 protected:
168 int ventricle_var_id; ///< Variable index of either left or right ventricle
169 int im_param_id; ///< Index of parameter Im
170};
171
172#endif // SVZERODSOLVER_MODEL_CLOSEDLOOPCORONARYBC_HPP_
model::Block source file
BlockType
The types of blocks supported by the solver.
Definition BlockType.h:15
BlockClass
The classes/categories of blocks supported. Some classes require special handling (e....
Definition BlockType.h:37
model::ClosedLoopHeartPulmonary source file
SparseSystem source file.
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.
Definition Block.h:100
const int id
Global ID of the block.
Definition Block.h:77
const Model * model
The model to which the block belongs.
Definition Block.h:78
const BlockType block_type
Type of this block.
Definition Block.h:79
TripletsContributions num_triplets
Number of triplets of element.
Definition ClosedLoopCoronaryBC.h:165
ParamId
Local IDs of the parameters.
Definition ClosedLoopCoronaryBC.h:118
void update_constant(SparseSystem &system, std::vector< double > &parameters)
Update the constant contributions of the element in a sparse system.
Definition ClosedLoopCoronaryBC.cpp:11
void setup_dofs(DOFHandler &dofhandler)
Set up the degrees of freedom (DOF) of the block.
Definition ClosedLoopCoronaryBC.cpp:7
void update_solution(SparseSystem &system, std::vector< double > &parameters, 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.
Definition ClosedLoopCoronaryBC.cpp:36
int ventricle_var_id
Variable index of either left or right ventricle.
Definition ClosedLoopCoronaryBC.h:168
ClosedLoopCoronaryBC(int id, Model *model, BlockType block_type)
Construct a ClosedLoopCoronaryBC object.
Definition ClosedLoopCoronaryBC.h:106
virtual void setup_model_dependent_params()=0
Setup parameters that depend on the model.
int im_param_id
Index of parameter Im.
Definition ClosedLoopCoronaryBC.h:169
Degree-of-freedom handler.
Definition DOFHandler.h:21
Model of 0D elements.
Definition Model.h:48
Sparse system.
Definition SparseSystem.h:30
Handles the properties of input parameters.
Definition Parameter.h:100
The number of triplets that the element contributes to the global system.
Definition Block.h:26