svZeroDSolver
Loading...
Searching...
No Matches
BloodVesselCRL.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/**
5 * @file BloodVesselCRL.h
6 * @brief model::BloodVesselCRL source file
7 */
8#ifndef SVZERODSOLVER_MODEL_BLOODVESSELCRL_HPP_
9#define SVZERODSOLVER_MODEL_BLOODVESSELCRL_HPP_
10
11#include <math.h>
12
13#include "Block.h"
14#include "SparseSystem.h"
15
16/**
17 * @brief Capacitor-resistor-inductor blood vessel with optional stenosis
18 *
19 * Models the mechanical behavior of a bloodvesselCRL with optional stenosis.
20 *
21 * \f[
22 * \begin{circuitikz} \draw
23 * node[left] {$Q_{in}$} [-latex] (0,0) -- (0.8,0);
24 * \draw (1,0) node[anchor=south]{$P_{in}$}
25 * to [R, l=$R$, *-] (3,0)
26 * to [R, l=$S$, -] (5,0)
27 * (5,0) to [L, l=$L$, -*] (7,0)
28 * node[anchor=south]{$P_{out}$}
29 * (1,0) to [C, l=$C$, -] (1,-1.5)
30 * node[ground]{};
31 * \draw [-latex] (7.2,0) -- (8,0) node[right] {$Q_{out}$};
32 * \end{circuitikz}
33 * \f]
34 *
35 * ### Governing equations
36 *
37 * \f[
38 * P_\text{in}-P_\text{out} - (R + S|Q_\text{out}|) Q_\text{out}-L
39 * \dot{Q}_\text{out}=0 \f]
40 *
41 * \f[
42 * Q_\text{in}-Q_\text{out} - C \dot{P}_\text{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 & 0 & -1 & -R \\
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 & 0 & 0 & 0
61 * \end{array}\right]
62 * \f]
63 *
64 * \f[
65 * \mathbf{c}^{e} = S|Q_\text{out}|
66 * \left[\begin{array}{c}
67 * -Q_\text{out} \\
68 * 0
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{out} & 0 & 0 \\
77 * 0 & 0 & 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{out}|
84 * \left[\begin{array}{cccc}
85 * 0 & 0 & 0 & 0 \\
86 * 0 & 0 & 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_3 & 0 & -\dot{y}_3 & -|y_3|y_3 \\
102 * 0 & 0 & -\dot{y}_0 & 0 \\
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 */
116class BloodVesselCRL : public Block {
117 public:
118 /**
119 * @brief Local IDs of the parameters
120 *
121 */
122 enum ParamId {
123 RESISTANCE = 0,
124 CAPACITANCE = 1,
125 INDUCTANCE = 2,
126 STENOSIS_COEFFICIENT = 3,
127 };
128
129 /**
130 * @brief Construct a new BloodVesselCRL object
131 *
132 * @param id Global ID of the block
133 * @param model The model to which the block belongs
134 */
136 : Block(id, model, BlockType::blood_vessel_CRL, BlockClass::vessel,
137 {{"R_poiseuille", InputParameter()},
138 {"C", InputParameter(true)},
139 {"L", InputParameter(true)},
140 {"stenosis_coefficient", InputParameter(true)}}) {}
141
142 /**
143 * @brief Set up the degrees of freedom (DOF) of the block
144 *
145 * Set \ref global_var_ids and \ref global_eqn_ids of the element based on the
146 * number of equations and the number of internal variables of the
147 * element.
148 *
149 * @param dofhandler Degree-of-freedom handler to register variables and
150 * equations at
151 */
152 void setup_dofs(DOFHandler& dofhandler);
153
154 /**
155 * @brief Update the constant contributions of the element in a sparse
156 system
157 *
158 * @param system System to update contributions at
159 * @param parameters Parameters of the model
160 */
161 void update_constant(SparseSystem& system, std::vector<double>& parameters);
162
163 /**
164 * @brief Update the solution-dependent contributions of the element in a
165 * sparse system
166 *
167 * @param system System to update contributions at
168 * @param parameters Parameters of the model
169 * @param y Current solution
170 * @param dy Current derivate of the solution
171 */
172 void update_solution(SparseSystem& system, std::vector<double>& parameters,
173 const Eigen::Matrix<double, Eigen::Dynamic, 1>& y,
174 const Eigen::Matrix<double, Eigen::Dynamic, 1>& dy);
175
176 /**
177 * @brief Set the gradient of the block contributions with respect to the
178 * parameters
179 *
180 * @param jacobian Jacobian with respect to the parameters
181 * @param alpha Current parameter vector
182 * @param residual Residual with respect to the parameters
183 * @param y Current solution
184 * @param dy Time-derivative of the current solution
185 */
186 void update_gradient(Eigen::SparseMatrix<double>& jacobian,
187 Eigen::Matrix<double, Eigen::Dynamic, 1>& residual,
188 Eigen::Matrix<double, Eigen::Dynamic, 1>& alpha,
189 std::vector<double>& y, std::vector<double>& dy);
190
191 /**
192 * @brief Number of triplets of element
193 *
194 * Number of triplets that the element contributes to the global system
195 * (relevant for sparse memory reservation)
196 */
198};
199
200#endif // SVZERODSOLVER_MODEL_BLOODVESSELCRL_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
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 BloodVesselCRL.cpp:10
ParamId
Local IDs of the parameters.
Definition BloodVesselCRL.h:122
BloodVesselCRL(int id, Model *model)
Construct a new BloodVesselCRL object.
Definition BloodVesselCRL.h:135
TripletsContributions num_triplets
Number of triplets of element.
Definition BloodVesselCRL.h:197
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 BloodVesselCRL.cpp:31
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 BloodVesselCRL.cpp:51
void setup_dofs(DOFHandler &dofhandler)
Set up the degrees of freedom (DOF) of the block.
Definition BloodVesselCRL.cpp:6
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