svZeroDSolver
Loading...
Searching...
No Matches
BloodVesselJunction.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 BloodVesselJunction.h
5 * @brief model::BloodVesselJunction source file
6 */
7#ifndef SVZERODSOLVER_MODEL_BLOODVESSELJUNCTION_HPP_
8#define SVZERODSOLVER_MODEL_BLOODVESSELJUNCTION_HPP_
9
10#include "Block.h"
11#include "BloodVessel.h"
12#include "BloodVesselCRL.h"
13#include "SparseSystem.h"
14
15/**
16 * @brief Junction between blood vessels
17 *
18 * Models a junction with one inlet and arbitrary outlets using
19 * modified blood vessel elements between each inlet and outlet pair.
20 *
21 * \f[
22 * \begin{circuitikz}
23 * \draw node[left] {$Q_\text{in}$} [-latex] (0,0) -- (0.8,0);
24 * \draw (1,0.1) node[anchor=south]{$P_\text{in}$};
25 * \draw (1,0) to [short, *-] (2.5,0.75);
26 * \draw (1,0) to [short, *-] (2.5,-0.75);
27 * \draw (2.5,0.75) node[anchor=south]{} to [generic, l_=$BV_{1}$, -*]
28 * (4.5,0.75); \draw (2.4,0.75) node[anchor=south]{}; \draw (4.6,0.75)
29 * node[anchor=south] {$P_{out,1}$}; \draw (2.5,-0.75) node[anchor=south]{}
30 to
31 * [generic, l^=$BV_{2}$, -*] (4.5,-0.75); \draw (2.4,-0.75)
32 * node[anchor=north]{}; \draw (4.6,-0.75) node[anchor=north]
33 * {$P_{out,2}$}; \draw [-latex] (4.7,0.75) -- (5.5,0.75) node[right]
34 * {$Q_{out,1}$}; \draw [-latex] (4.7,-0.75) -- (5.5,-0.75) node[right]
35 * {$Q_{out,2}$}; \end{circuitikz}
36 * \f]
37 *
38 * Each blood vessel is modelled as:
39 *
40 * \f[
41 * \begin{circuitikz} \draw
42 * node[left] {$Q_\text{in}$} [-latex] (0,0) -- (0.8,0);
43 * \draw (1,0) node[anchor=south]{$P_\text{in}$}
44 * to [R, l=$R$, *-] (3,0)
45 * to [R, l=$S$, -] (5,0)
46 * (5,0) to [L, l=$L$, -*] (7,0)
47 * node[anchor=south]{$P_\text{out}$};
48 * \draw [-latex] (7.2,0) -- (8,0) node[right] {$Q_\text{out}$};
49 * \end{circuitikz}
50 * \f]
51 *
52 * ### Governing equations
53 *
54 * \f[
55 * Q_\text{in}-\sum_{i}^{n_{outlets}} Q_{out, i} = 0
56 * \f]
57 *
58 * \f[
59 * P_\text{in}-P_{\text{out},i} - (R+S|Q_{\text{out},i}|) \cdot Q_{\text{out},i}
60 * - L \dot{Q}_{\text{out},i} = 0 \quad \forall i \in n_{outlets} \f]
61 *
62 * ### Local contributions
63 *
64 * \f[
65 * \mathbf{y}^{e}=\left[\begin{array}{lllllll}P_\text{in} & Q_\text{in}
66 * & P_{out, 1} & Q_{out, 1} &
67 * \dots & P_{out, i} & Q_{out, i}\end{array}\right] \f]
68 *
69 * \f[
70 * \mathbf{F}^{e} = \left[\begin{array}{ccccccc}
71 * 0 & 1 & 0 & -1 & 0 & -1 & \dots \\
72 * 1 & 0 & -1 & -R_1 & 0 & 0 & \\
73 * 1 & 0 & 0 & 0 & -1 & -R_2 & \\
74 * \vdots & & & & & & \ddots
75 * \end{array}\right]
76 * \f]
77 *
78 * \f[
79 * \mathbf{E}^{e} = \left[\begin{array}{ccccccc}
80 * 0 & 0 & 0 & 0 & 0 & 0 & \dots \\
81 * 0 & 0 & 0 & -L_1 & 0 & 0 & \\
82 * 0 & 0 & 0 & 0 & 0 & -L_2 & \\
83 * \vdots & & & & & & \ddots
84 * \end{array}\right]
85 * \f]
86 *
87 * \f[
88 * \mathbf{c}^{e} = \left[\begin{array}{c}
89 * 0 \\
90 * - S_1 |Q_{\text{in},1}| Q_{\text{in},1} \\
91 * - S_2 |Q_{\text{in},2}| Q_{\text{in},2} \\
92 * \vdots
93 * \end{array}\right] \f]
94 *
95 * ### Gradient
96 *
97 * Gradient of the equations with respect to the parameters:
98 *
99 * \f[
100 * \mathbf{J}^{e} = \left[\begin{array}{lllllllll}
101 * 0 & 0 & \dots & 0 & 0 & \dots & 0 & 0 &
102 * \dots \\
103 * - y_4 & 0 & \dots & - \dot y_4 & 0 & \dots & |y_4| y_4 & 0 & \dots \\
104 * 0 & - y_6 & \dots & 0 & - \dot y_6 & \dots & 0 & |y_6| y_6 & \dots \\
105 * 0 & 0 & \ddots & 0 & 0 & \ddots & 0 & 0 & \ddots \\
106 * \end{array}\right]
107 * \f]
108 *
109 * ### Parameters
110 *
111 * Parameter sequence for constructing this block
112 *
113 * * `i` Poiseuille resistance for inner blood vessel `i`
114 * * `i+num_outlets` Inductance for inner blood vessel `i`
115 * * `i+2*num_outlets` Stenosis coefficient for inner blood vessel `i`
116 *
117 * ### Usage in json configuration file
118 *
119 * "junctions": [
120 * {
121 * "inlet_vessels": [
122 * 0
123 * ],
124 * "junction_name": "J0",
125 * "junction_type": "BloodVesselJunction",
126 * "junction_values": {
127 * "R_poiseuille": [
128 * 100,
129 * 200
130 * ],
131 * "C": [
132 * 0.0,
133 * 0.0
134 * ],
135 * "L": [
136 * 0.0,
137 * 0.0
138 * ],
139 * "stenosis_coefficient": [
140 * 0.0,
141 * 0.0
142 * ]
143 * },
144 * "outlet_vessels": [
145 * 1,
146 * 2
147 * ]
148 * }
149 * ]
150 *
151 * ### Internal variables
152 *
153 * This block has no internal variables.
154 *
155 */
157 public:
158 /**
159 * @brief Construct a new BloodVesselJunction object
160 *
161 * @param id Global ID of the block
162 * @param model The model to which the block belongs
163 */
165 : Block(id, model, BlockType::blood_vessel_junction, BlockClass::junction,
166 {{"R_poiseuille", InputParameter()},
167 {"L", InputParameter()},
168 {"stenosis_coefficient", InputParameter()}}) {
169 input_params_list = true;
170 }
171
172 /**
173 * @brief Set up the degrees of freedom (DOF) of the block
174 *
175 * Set \ref global_var_ids and \ref global_eqn_ids of the element based on the
176 * number of equations and the number of internal variables of the
177 * element.
178 *
179 * @param dofhandler Degree-of-freedom handler to register variables and
180 * equations at
181 */
182 void setup_dofs(DOFHandler& dofhandler);
183
184 /**
185 * @brief Update the constant contributions of the element in a sparse system
186 *
187 * @param system System to update contributions at
188 * @param parameters Parameters of the model
189 */
190 void update_constant(SparseSystem& system, std::vector<double>& parameters);
191
192 /**
193 * @brief Update the solution-dependent contributions of the element in a
194 * sparse system
195 *
196 * @param system System to update contributions at
197 * @param parameters Parameters of the model
198 * @param y Current solution
199 * @param dy Current derivate of the solution
200 */
201 virtual void update_solution(
202 SparseSystem& system, std::vector<double>& parameters,
203 const Eigen::Matrix<double, Eigen::Dynamic, 1>& y,
204 const Eigen::Matrix<double, Eigen::Dynamic, 1>& dy);
205
206 /**
207 * @brief Set the gradient of the block contributions with respect to the
208 * parameters
209 *
210 * @param jacobian Jacobian with respect to the parameters
211 * @param alpha Current parameter vector
212 * @param residual Residual with respect to the parameters
213 * @param y Current solution
214 * @param dy Time-derivative of the current solution
215 */
216 void update_gradient(Eigen::SparseMatrix<double>& jacobian,
217 Eigen::Matrix<double, Eigen::Dynamic, 1>& residual,
218 Eigen::Matrix<double, Eigen::Dynamic, 1>& alpha,
219 std::vector<double>& y, std::vector<double>& dy);
220
221 /**
222 * @brief Number of triplets of element
223 *
224 * Number of triplets that the element contributes to the global system
225 * (relevant for sparse memory reservation)
226 */
228
229 private:
230 int num_outlets;
231};
232
233#endif // SVZERODSOLVER_MODEL_BLOODVESSELJUNCTION_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:41
model::BloodVessel source file
model::BloodVesselCRL source file
SparseSystem source file.
bool input_params_list
Are input parameters given as a list?
Definition Block.h:89
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
void setup_dofs(DOFHandler &dofhandler)
Set up the degrees of freedom (DOF) of the block.
Definition BloodVesselJunction.cpp:6
void update_constant(SparseSystem &system, std::vector< double > &parameters)
Update the constant contributions of the element in a sparse system.
Definition BloodVesselJunction.cpp:19
TripletsContributions num_triplets
Number of triplets of element.
Definition BloodVesselJunction.h:227
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.
Definition BloodVesselJunction.cpp:55
BloodVesselJunction(int id, Model *model)
Construct a new BloodVesselJunction object.
Definition BloodVesselJunction.h:164
virtual 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 BloodVesselJunction.cpp:38
Model of 0D elements.
Definition Model.h:52
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