svZeroDSolver
Loading...
Searching...
No Matches
ClosedLoopHeartPulmonary.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 ClosedLoopHeartPulmonary.h
5 * @brief model::ClosedLoopHeartPulmonary source file
6 */
7#ifndef SVZERODSOLVER_MODEL_CLOSEDLOOPHEARTPULMONARY_HPP_
8#define SVZERODSOLVER_MODEL_CLOSEDLOOPHEARTPULMONARY_HPP_
9
10#include "Block.h"
11#include "SparseSystem.h"
12
13/**
14 * @brief Define math constants (for M_PI)
15 *
16 */
17#define _USE_MATH_DEFINES
18
19#include <cmath>
20
21/**
22 * @brief Heart and pulmonary circulation model
23 *
24 * Models the mechanics of the 4 heart chambers and pulmonary circulation
25 *
26 * References: \cite sankaran2012patient and \cite menon2023predictors
27 *
28 * TODO: Equations and circuit diagram
29 *
30 * ### Parameters
31 *
32 * Parameter sequence for constructing this block
33 *
34 * * `0` Atrial systole time fraction
35 * * `1` Time for P-wave
36 * * `2` Scaling for right ventricle elastance
37 * * `3` Scaling for left ventricle elastance
38 * * `4` Scaling for intramyocardial pressure (left coronaries)
39 * * `5` Scaling for intramyocardial pressure (right coronaries)
40 * * `6` Right atrium inductance
41 * * `7` Right atrium outflow resistance
42 * * `8` Right ventricle inductance
43 * * `9` Right ventricle outflow resistance
44 * * `10` Left atrium inductance
45 * * `11` Left atrium outflow resistance
46 * * `12` Left ventricle inductance
47 * * `13` Left ventricle outflow resistance
48 * * `14` Right ventricle unstressed volume
49 * * `15` Left ventricle unstressed volume
50 * * `16` Pulmonary resistance
51 * * `17` Pulmonary capacitance
52 * * `18` Aortic capacitance
53 * * `19` Right atrium pressure scaling
54 * * `20` Right atrium volume scaling
55 * * `21` Left atrium pressure scaling
56 * * `22` Left atrium volume scaling
57 * * `23` Right atrium elastance
58 * * `24` Left atrium elastance
59 * * `25` Right atrium resting volume
60 * * `26` Left atrium resting volume
61 *
62 * ### Usage in json configuration file
63 *
64 * "closed_loop_blocks": [
65 * {
66 * "outlet_blocks": [
67 * "branch0_seg0"
68 * ],
69 * "closed_loop_type": "ClosedLoopHeartAndPulmonary",
70 * "cardiac_cycle_period": 1.0169,
71 * "parameters": {
72 * "Tsa": 0.40742,
73 * "tpwave": 8.976868,
74 * "Erv_s": 2.125279,
75 * "Elv_s": 3.125202,
76 * "iml": 0.509365,
77 * "imr": 0.806369,
78 * "Lrv_a": 0.000186865,
79 * "Rrv_a": 0.035061704,
80 * "Lra_v": 0.000217032,
81 * "Rra_v": 0.007887459,
82 * "Lla_v": 0.000351787,
83 * "Rla_v": 0.005310825,
84 * "Rlv_ao": 0.034320234,
85 * "Llv_a": 0.000110776,
86 * "Vrv_u": 9.424629,
87 * "Vlv_u": 5.606007,
88 * "Rpd": 0.098865401,
89 * "Cp": 1.090989,
90 * "Cpa": 0.556854,
91 * "Kxp_ra": 9.22244,
92 * "Kxv_ra": 0.004837,
93 * "Emax_ra": 0.208858,
94 * "Vaso_ra": 4.848742,
95 * "Kxp_la": 9.194992,
96 * "Kxv_la": 0.008067,
97 * "Emax_la": 0.303119,
98 * "Vaso_la": 9.355754
99 * }
100 * }
101 * ],
102 * "initial_condition": {
103 * "V_RA:CLH": 38.43,
104 * "V_RV:CLH": 96.07,
105 * "V_LA:CLH": 38.43,
106 * "V_LV:CLH": 96.07,
107 * "P_pul:CLH": 8.0
108 * }
109 *
110 * ### Internal variables
111 *
112 * Names of internal variables in this block's output:
113 *
114 * * `V_RA`: Right atrium volume
115 * * `Q_RA`: Right atrium outflow
116 * * `P_RV`: Right ventricle pressure
117 * * `V_RV`: Right ventricle volume
118 * * `Q_RV`: Right ventricle outflow
119 * * `P_pul`: Pulmonary pressure
120 * * `P_LA`: Left atrium pressure
121 * * `V_LA`: Left atrium volume
122 * * `Q_LA`: Left atrium outflow
123 * * `P_LV`: Left ventricle pressure
124 * * `V_LV`: Left ventricle volume
125 * * `Q_LV`: Left ventricle outflow
126 *
127 */
129 public:
130 /**
131 * @brief Construct a new ClosedLoopHeartPulmonary object
132 *
133 * @param id Global ID of the block
134 * @param model The model to which the block belongs
135 */
137 : Block(id, model, BlockType::closed_loop_heart_pulmonary,
138 BlockClass::closed_loop,
139 {{"Tsa", InputParameter()}, {"tpwave", InputParameter()},
140 {"Erv_s", InputParameter()}, {"Elv_s", InputParameter()},
141 {"iml", InputParameter()}, {"imr", InputParameter()},
142 {"Lra_v", InputParameter()}, {"Rra_v", InputParameter()},
143 {"Lrv_a", InputParameter()}, {"Rrv_a", InputParameter()},
144 {"Lla_v", InputParameter()}, {"Rla_v", InputParameter()},
145 {"Llv_a", InputParameter()}, {"Rlv_ao", InputParameter()},
146 {"Vrv_u", InputParameter()}, {"Vlv_u", InputParameter()},
147 {"Rpd", InputParameter()}, {"Cp", InputParameter()},
148 {"Cpa", InputParameter()}, {"Kxp_ra", InputParameter()},
149 {"Kxv_ra", InputParameter()}, {"Kxp_la", InputParameter()},
150 {"Kxv_la", InputParameter()}, {"Emax_ra", InputParameter()},
151 {"Emax_la", InputParameter()}, {"Vaso_ra", InputParameter()},
152 {"Vaso_la", InputParameter()}}) {}
153
154 /**
155 * @brief Local IDs of the parameters
156 *
157 */
158 enum ParamId {
159 TSA = 0, ///< Fractions of cardiac cycle (not sure)
160 TPWAVE = 1, ///< Fraction of cardiac cycle (P-wave)
161 ERV_S = 2, ///< Scaling for right ventricle elastance
162 ELV_S = 3, ///< Scaling for left ventricle elastance
163 IML = 4, ///< Scaling for intramyocardial pressure (left coronaries)
164 IMR = 5, ///< Scaling for intramyocardial pressure (right coronaries)
165 LRA_V = 6, ///< Right atrium inductance
166 RRA_V = 7, ///< Right atrium outflow resistance
167 LRV_A = 8, ///< Right ventricle inductance
168 RRV_A = 9, ///< Right ventricle outflow resistance
169 LLA_V = 10, ///< Left atrium inductance
170 RLA_V = 11, ///< Left atrium outflow resistance
171 LLV_A = 12, ///< Left ventricle inductance
172 RLV_AO = 13, ///< Left ventricle outflow resistance
173 VRV_U = 14, ///< Right ventricle unstressed volume
174 VLV_U = 15, ///< Left ventricle unstressed volume
175 RPD = 16, ///< Pulmonary resistance
176 CP = 17, ///< Pulmonary capacitance
177 CPA = 18, ///< Aortic capacitance
178 KXP_RA = 19, ///< Right atrium pressure-volume relationship (?)
179 KXV_RA = 20, ///< Right atrium pressure-volume relationship (?)
180 KXP_LA = 21, ///< Left atrium pressure-volume relationship (?)
181 KXV_LA = 22, ///< Left atrium pressure-volume relationship (?)
182 EMAX_RA = 23, ///< Right atrium elastance (?)
183 EMAX_LA = 24, ///< Left atrium elastance (?)
184 VASO_RA = 25, ///< Right atrium rest volume (?)
185 VASO_LA = 26, ///< Left atrium rest volume (?)
186 };
187
188 /**
189 * @brief Set up the degrees of freedom (DOF) of the block
190 *
191 * Set \ref global_var_ids and \ref global_eqn_ids of the element based on the
192 * number of equations and the number of internal variables of the
193 * element.
194 *
195 * @param dofhandler Degree-of-freedom handler to register variables and
196 * equations at
197 */
198 void setup_dofs(DOFHandler& dofhandler);
199
200 /**
201 * @brief Update the constant contributions of the element in a sparse
202 system
203 *
204 * @param system System to update contributions at
205 * @param parameters Parameters of the model
206 */
207 void update_constant(SparseSystem& system, std::vector<double>& parameters);
208
209 /**
210 * @brief Update the time-dependent contributions of the element in a sparse
211 * system
212 *
213 * @param system System to update contributions at
214 * @param parameters Parameters of the model
215 */
216 void update_time(SparseSystem& system, std::vector<double>& parameters);
217
218 /**
219 * @brief Update the solution-dependent contributions of the element in a
220 * sparse system
221 *
222 * @param system System to update contributions at
223 * @param parameters Parameters of the model
224 * @param y Current solution
225 * @param dy Current derivate of the solution
226 */
227 void update_solution(SparseSystem& system, std::vector<double>& parameters,
228 const Eigen::Matrix<double, Eigen::Dynamic, 1>& y,
229 const Eigen::Matrix<double, Eigen::Dynamic, 1>& dy);
230
231 /**
232 * @brief Modify the solution after solving it
233 *
234 * @param y Current solution
235 */
236 void post_solve(Eigen::Matrix<double, Eigen::Dynamic, 1>& y);
237
238 /**
239 * @brief Number of triplets of element
240 *
241 * Number of triplets that the element contributes to the global system
242 * (relevant for sparse memory reservation)
243 */
245
246 private:
247 // Below variables change every timestep and are then combined with
248 // expressions that are updated with solution
249 double AA = 0.0; // Atrial activation function
250 double Elv = 0.0; // LV elastance
251 double Erv = 0.0; // RV elastance
252 double psi_ra, psi_la, psi_ra_derivative,
253 psi_la_derivative; // Expressions for atrial activation
254 double valves[16];
255
256 /**
257 * @brief Update the atrial activation and LV/RV elastance functions which
258 * depend on time
259 *
260 * @param parameters Parameters of the model
261 */
262 void get_activation_and_elastance_functions(std::vector<double>& parameters);
263
264 /**
265 * @brief Compute sub-expressions that are part of atrial elastance and
266 * depends on atrial volume from the solution vector
267 *
268 * @param parameters Parameters of the model
269 * @param y Current solution
270 */
271 void get_psi_ra_la(std::vector<double>& parameters,
272 const Eigen::Matrix<double, Eigen::Dynamic, 1>& y);
273
274 /**
275 * @brief Valve positions for each heart chamber
276 *
277 * @param y Current solution
278 */
279 void get_valve_positions(const Eigen::Matrix<double, Eigen::Dynamic, 1>& y);
280};
281
282#endif // SVZERODSOLVER_MODEL_CLOSEDLOOPHEARTPULMONARY_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:38
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
ClosedLoopHeartPulmonary(int id, Model *model)
Construct a new ClosedLoopHeartPulmonary object.
Definition ClosedLoopHeartPulmonary.h:136
void setup_dofs(DOFHandler &dofhandler)
Set up the degrees of freedom (DOF) of the block.
Definition ClosedLoopHeartPulmonary.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 ClosedLoopHeartPulmonary.cpp:101
void post_solve(Eigen::Matrix< double, Eigen::Dynamic, 1 > &y)
Modify the solution after solving it.
Definition ClosedLoopHeartPulmonary.cpp:282
ParamId
Local IDs of the parameters.
Definition ClosedLoopHeartPulmonary.h:158
@ KXV_RA
Right atrium pressure-volume relationship (?)
Definition ClosedLoopHeartPulmonary.h:179
@ RRA_V
Right atrium outflow resistance.
Definition ClosedLoopHeartPulmonary.h:166
@ KXP_RA
Right atrium pressure-volume relationship (?)
Definition ClosedLoopHeartPulmonary.h:178
@ VRV_U
Right ventricle unstressed volume.
Definition ClosedLoopHeartPulmonary.h:173
@ LRV_A
Right ventricle inductance.
Definition ClosedLoopHeartPulmonary.h:167
@ RLV_AO
Left ventricle outflow resistance.
Definition ClosedLoopHeartPulmonary.h:172
@ EMAX_RA
Right atrium elastance (?)
Definition ClosedLoopHeartPulmonary.h:182
@ TPWAVE
Fraction of cardiac cycle (P-wave)
Definition ClosedLoopHeartPulmonary.h:160
@ IML
Scaling for intramyocardial pressure (left coronaries)
Definition ClosedLoopHeartPulmonary.h:163
@ IMR
Scaling for intramyocardial pressure (right coronaries)
Definition ClosedLoopHeartPulmonary.h:164
@ VASO_LA
Left atrium rest volume (?)
Definition ClosedLoopHeartPulmonary.h:185
@ LLV_A
Left ventricle inductance.
Definition ClosedLoopHeartPulmonary.h:171
@ ELV_S
Scaling for left ventricle elastance.
Definition ClosedLoopHeartPulmonary.h:162
@ CP
Pulmonary capacitance.
Definition ClosedLoopHeartPulmonary.h:176
@ VLV_U
Left ventricle unstressed volume.
Definition ClosedLoopHeartPulmonary.h:174
@ RLA_V
Left atrium outflow resistance.
Definition ClosedLoopHeartPulmonary.h:170
@ EMAX_LA
Left atrium elastance (?)
Definition ClosedLoopHeartPulmonary.h:183
@ TSA
Fractions of cardiac cycle (not sure)
Definition ClosedLoopHeartPulmonary.h:159
@ KXP_LA
Left atrium pressure-volume relationship (?)
Definition ClosedLoopHeartPulmonary.h:180
@ CPA
Aortic capacitance.
Definition ClosedLoopHeartPulmonary.h:177
@ RPD
Pulmonary resistance.
Definition ClosedLoopHeartPulmonary.h:175
@ LLA_V
Left atrium inductance.
Definition ClosedLoopHeartPulmonary.h:169
@ ERV_S
Scaling for right ventricle elastance.
Definition ClosedLoopHeartPulmonary.h:161
@ VASO_RA
Right atrium rest volume (?)
Definition ClosedLoopHeartPulmonary.h:184
@ KXV_LA
Left atrium pressure-volume relationship (?)
Definition ClosedLoopHeartPulmonary.h:181
@ RRV_A
Right ventricle outflow resistance.
Definition ClosedLoopHeartPulmonary.h:168
@ LRA_V
Right atrium inductance.
Definition ClosedLoopHeartPulmonary.h:165
void update_time(SparseSystem &system, std::vector< double > &parameters)
Update the time-dependent contributions of the element in a sparse system.
Definition ClosedLoopHeartPulmonary.cpp:78
TripletsContributions num_triplets
Number of triplets of element.
Definition ClosedLoopHeartPulmonary.h:244
void update_constant(SparseSystem &system, std::vector< double > &parameters)
Update the constant contributions of the element in a sparse system.
Definition ClosedLoopHeartPulmonary.cpp:13
Degree-of-freedom handler.
Definition DOFHandler.h:21
Model of 0D elements.
Definition Model.h:49
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