svZeroDSolver
Loading...
Searching...
No Matches
BloodVessel.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 BloodVessel.h
5 * @brief model::BloodVessel source file
6 */
7#ifndef SVZERODSOLVER_MODEL_BLOODVESSEL_HPP_
8#define SVZERODSOLVER_MODEL_BLOODVESSEL_HPP_
9
10#include <math.h>
11
12#include "Block.h"
13#include "SparseSystem.h"
14
15/**
16 * @brief Resistor-capacitor-inductor blood vessel with optional stenosis
17 *
18 * Models the mechanical behavior of a bloodvessel with optional stenosis.
19 *
20 * \f[
21 * \begin{circuitikz} \draw
22 * node[left] {$Q_{in}$} [-latex] (0,0) -- (0.8,0);
23 * \draw (1,0) node[anchor=south]{$P_{in}$}
24 * to [R, l=$R$, *-] (3,0)
25 * to [R, l=$S$, -] (5,0)
26 * (5,0) to [L, l=$L$, -*] (7,0)
27 * node[anchor=south]{$P_{out}$}
28 * (5,0) to [C, l=$C$, -] (5,-1.5)
29 * node[ground]{};
30 * \draw [-latex] (7.2,0) -- (8,0) node[right] {$Q_{out}$};
31 * \end{circuitikz}
32 * \f]
33 *
34 * ### Governing equations
35 *
36 * \f[
37 * P_\text{in}-P_\text{out} - (R + S|Q_\text{in}|) Q_\text{in}-L
38 * \dot{Q}_\text{out}=0 \f]
39 *
40 * \f[
41 * Q_\text{in}-Q_\text{out} - C \dot{P}_\text{in}+C(R +
42 * 2S|Q_\text{in}|) \dot{Q}_{in}=0 \f]
43 *
44 * ### Local contributions
45 *
46 * \f[
47 * \mathbf{y}^{e}=\left[\begin{array}{llll}P_{i n} & Q_{in} &
48 * P_{out} & Q_{out}\end{array}\right]^\text{T} \f]
49 *
50 * \f[
51 * \mathbf{F}^{e}=\left[\begin{array}{cccc}
52 * 1 & -R & -1 & 0 \\
53 * 0 & 1 & 0 & -1
54 * \end{array}\right]
55 * \f]
56 *
57 * \f[
58 * \mathbf{E}^{e}=\left[\begin{array}{cccc}
59 * 0 & 0 & 0 & -L \\
60 * -C & CR & 0 & 0
61 * \end{array}\right]
62 * \f]
63 *
64 * \f[
65 * \mathbf{c}^{e} = S|Q_\text{in}|
66 * \left[\begin{array}{c}
67 * -Q_\text{in} \\
68 * 2C\dot{Q}_\text{in}
69 * \end{array}\right]
70 * \f]
71 *
72 * \f[
73 * \left(\frac{\partial\mathbf{c}}{\partial\mathbf{y}}\right)^{e} =
74 * S \text{sgn} (Q_\text{in})
75 * \left[\begin{array}{cccc}
76 * 0 & -2Q_\text{in} & 0 & 0 \\
77 * 0 & 2C\dot{Q}_\text{in} & 0 & 0
78 * \end{array}\right]
79 * \f]
80 *
81 * \f[
82 * \left(\frac{\partial\mathbf{c}}{\partial\dot{\mathbf{y}}}\right)^{e} =
83 * S|Q_\text{in}|
84 * \left[\begin{array}{cccc}
85 * 0 & 0 & 0 & 0 \\
86 * 0 & 2C & 0 & 0
87 * \end{array}\right]
88 * \f]
89 *
90 * with the stenosis resistance \f$ S=K_{t} \frac{\rho}{2
91 * A_{o}^{2}}\left(\frac{A_{o}}{A_{s}}-1\right)^{2} \f$.
92 * \f$R\f$, \f$C\f$, and \f$L\f$ refer to
93 * Poisieuille resistance, capacitance and inductance, respectively.
94 *
95 * ### Gradient
96 *
97 * Gradient of the equations with respect to the parameters:
98 *
99 * \f[
100 * \mathbf{J}^{e} = \left[\begin{array}{cccc}
101 * -y_2 & 0 & -\dot{y}_4 & -|y_2|y_2 \\
102 * C\dot{y}_2 & (-\dot{y}_1+(R+2S|Q_\text{in}|)\dot{y}_2) & 0 & 2C|y_2|\dot{y}_2
103 * \end{array}\right]
104 * \f]
105 *
106 * ### Parameters
107 *
108 * Parameter sequence for constructing this block
109 *
110 * * `0` Poiseuille resistance
111 * * `1` Capacitance
112 * * `2` Inductance
113 * * `3` Stenosis coefficient
114 *
115 * ### Internal variables
116 *
117 * This block has no internal variables.
118 *
119 */
120class BloodVessel : public Block {
121 public:
122 /**
123 * @brief Local IDs of the parameters
124 *
125 */
126 enum ParamId {
127 RESISTANCE = 0,
128 CAPACITANCE = 1,
129 INDUCTANCE = 2,
130 STENOSIS_COEFFICIENT = 3,
131 };
132
133 /**
134 * @brief Construct a new BloodVessel object
135 *
136 * @param id Global ID of the block
137 * @param model The model to which the block belongs
138 */
140 : Block(id, model, BlockType::blood_vessel, BlockClass::vessel,
141 {{"R_poiseuille", InputParameter()},
142 {"C", InputParameter(true)},
143 {"L", InputParameter(true)},
144 {"stenosis_coefficient", InputParameter(true)}}) {}
145
146 /**
147 * @brief Set up the degrees of freedom (DOF) of the block
148 *
149 * Set \ref global_var_ids and \ref global_eqn_ids of the element based on the
150 * number of equations and the number of internal variables of the
151 * element.
152 *
153 * @param dofhandler Degree-of-freedom handler to register variables and
154 * equations at
155 */
156 void setup_dofs(DOFHandler &dofhandler);
157
158 /**
159 * @brief Update the constant contributions of the element in a sparse
160 system
161 *
162 * @param system System to update contributions at
163 * @param parameters Parameters of the model
164 */
165 void update_constant(SparseSystem &system, std::vector<double> &parameters);
166
167 /**
168 * @brief Update the solution-dependent contributions of the element in a
169 * sparse system
170 *
171 * @param system System to update contributions at
172 * @param parameters Parameters of the model
173 * @param y Current solution
174 * @param dy Current derivate of the solution
175 */
176 void update_solution(SparseSystem &system, std::vector<double> &parameters,
177 const Eigen::Matrix<double, Eigen::Dynamic, 1> &y,
178 const Eigen::Matrix<double, Eigen::Dynamic, 1> &dy);
179
180 /**
181 * @brief Set the gradient of the block contributions with respect to the
182 * parameters
183 *
184 * @param jacobian Jacobian with respect to the parameters
185 * @param alpha Current parameter vector
186 * @param residual Residual with respect to the parameters
187 * @param y Current solution
188 * @param dy Time-derivative of the current solution
189 */
190 void update_gradient(Eigen::SparseMatrix<double> &jacobian,
191 Eigen::Matrix<double, Eigen::Dynamic, 1> &residual,
192 Eigen::Matrix<double, Eigen::Dynamic, 1> &alpha,
193 std::vector<double> &y, std::vector<double> &dy);
194
195 /**
196 * @brief Number of triplets of element
197 *
198 * Number of triplets that the element contributes to the global system
199 * (relevant for sparse memory reservation)
200 */
202};
203
204#endif // SVZERODSOLVER_MODEL_BLOODVESSEL_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
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
void update_constant(SparseSystem &system, std::vector< double > &parameters)
Update the constant contributions of the element in a sparse system.
Definition BloodVessel.cpp:10
ParamId
Local IDs of the parameters.
Definition BloodVessel.h:126
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 BloodVessel.cpp:33
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 BloodVessel.cpp:59
BloodVessel(int id, Model *model)
Construct a new BloodVessel object.
Definition BloodVessel.h:139
TripletsContributions num_triplets
Number of triplets of element.
Definition BloodVessel.h:201
void setup_dofs(DOFHandler &dofhandler)
Set up the degrees of freedom (DOF) of the block.
Definition BloodVessel.cpp:6
Model of 0D elements.
Definition Model.h:48
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