svZeroDSolver
Loading...
Searching...
No Matches
Solver.h
Go to the documentation of this file.
1// Copyright (c) Stanford University, The Regents of the University of
2// California, and others.
3//
4// All Rights Reserved.
5//
6// See Copyright-SimVascular.txt for additional details.
7//
8// Permission is hereby granted, free of charge, to any person obtaining
9// a copy of this software and associated documentation files (the
10// "Software"), to deal in the Software without restriction, including
11// without limitation the rights to use, copy, modify, merge, publish,
12// distribute, sublicense, and/or sell copies of the Software, and to
13// permit persons to whom the Software is furnished to do so, subject
14// to the following conditions:
15//
16// The above copyright notice and this permission notice shall be included
17// in all copies or substantial portions of the Software.
18//
19// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
20// IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
21// TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
22// PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER
23// OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
24// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
25// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
26// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
27// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
28// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
29// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
30/**
31 * @file Solver.h
32 * @brief Solver source file
33 */
34
35#include "Integrator.h"
36#include "Model.h"
38#include "State.h"
39#include "debug.h"
40
41#ifndef SVZERODSOLVER_SOLVE_SOLVER_HPP_
42#define SVZERODSOLVER_SOLVE_SOLVER_HPP_
43
44/**
45 * @brief Class for running 0D simulations.
46 *
47 * The solver solves for pressure and flow rate at the nodes of the
48 * lumped-parameter system \cite pfaller22. The lumped-parameter network is
49 * documented in SparseSystem. Refer to Integrator for the time discretization.
50 *
51 */
52class Solver {
53 public:
54 /**
55 * @brief Construct a new Solver object
56 *
57 * @param config Configuration handler
58 */
59 Solver(const nlohmann::json& config);
60
61 /**
62 * @brief Run the simulation
63 *
64 */
65 void run();
66
67 /**
68 * @brief Get the full result as a csv encoded string
69 *
70 * @return std::string Result
71 */
72 std::string get_full_result() const;
73
74 /**
75 * @brief Get the result of a single DOF over time
76 *
77 * @param dof_name Name of the degree-of-freedom
78 * @return Eigen::VectorXd Result
79 */
80 Eigen::VectorXd get_single_result(const std::string& dof_name) const;
81
82 /**
83 * @brief Get the result of a single DOF averaged over time
84 *
85 * @param dof_name Name of the degree-of-freedom
86 * @return double Result
87 */
88 double get_single_result_avg(const std::string& dof_name) const;
89
90 /**
91 * @brief Get the time steps of the result
92 *
93 * @return std::vector<double> Vector of times
94 */
95 std::vector<double> get_times() const;
96
97 /**
98 * @brief Update the parameters of a block
99 *
100 * @param block_name Name of the block
101 * @param new_params New parameters
102 */
103 void update_block_params(const std::string& block_name,
104 const std::vector<double>& new_params);
105
106 /**
107 * @brief Read the parameters of a block
108 *
109 * @param block_name Name of the block
110 *
111 * @return std::vector<double> Block parameters
112 */
113 std::vector<double> read_block_params(const std::string& block_name);
114
115 /**
116 * @brief Write the result to a csv file.
117 *
118 * @param filename
119 */
120 void write_result_to_csv(const std::string& filename) const;
121
122 private:
123 std::shared_ptr<Model> model;
124 SimulationParameters simparams;
125 std::vector<State> states;
126 std::vector<double> times;
127 State initial_state;
128
129 void sanity_checks();
130
131 /**
132 * @brief Get indices of flow and pressure degrees-of-freedom in solution
133 * vector for all vessel caps
134 *
135 * @return std::vector<std::pair<int, int>> Indices of flow and pressure
136 * degrees-of-freedom in solution vector for all vessel caps
137 */
138 std::vector<std::pair<int, int>> get_vessel_caps_dof_indices();
139
140 /**
141 * @brief Check if flows and pressures for all vessel caps have converged,
142 * based on cycle-to-cycle error for last two simulated cardiac cycles
143 *
144 * @param states_last_two_cycles Vector of solution states for last two
145 * simulated cardiac cycles
146 * @param vessel_caps_dof_indices Indices of flow and pressure
147 * degrees-of-freedom in solution vector for all vessel caps
148 *
149 * @return bool True if flows and pressures for all vessel caps have converged
150 */
151 bool check_vessel_cap_convergence(
152 const std::vector<State>& states_last_two_cycles,
153 const std::vector<std::pair<int, int>>& vessel_caps_dof_indices);
154
155 /**
156 * @brief Get cycle-to-cycle errors for flow and pressure for a single vessel
157 * cap
158 *
159 * @param states_last_two_cycles Vector of solution states for last two
160 * simulated cardiac cycles
161 * @param dof_indices Indices of flow and pressure degrees-of-freedom in
162 * solution vector for a single vessel cap
163 *
164 * @return std::pair<double, double> Cycle-to-cycle errors for flow and
165 * pressure
166 */
167 std::pair<double, double> get_cycle_to_cycle_errors_in_flow_and_pressure(
168 const std::vector<State>& states_last_two_cycles,
169 const std::pair<int, int>& dof_indices);
170};
171
172#endif
Integrator source file.
model::Model source file
Source file to read simulation configuration.
State source file.
Class for running 0D simulations.
Definition Solver.h:52
std::vector< double > read_block_params(const std::string &block_name)
Read the parameters of a block.
Definition Solver.cpp:364
void run()
Run the simulation.
Definition Solver.cpp:40
Solver(const nlohmann::json &config)
Construct a new Solver object.
Definition Solver.cpp:5
void update_block_params(const std::string &block_name, const std::vector< double > &new_params)
Update the parameters of a block.
Definition Solver.cpp:340
double get_single_result_avg(const std::string &dof_name) const
Get the result of a single DOF averaged over time.
Definition Solver.cpp:328
std::vector< double > get_times() const
Get the time steps of the result.
Definition Solver.cpp:297
Eigen::VectorXd get_single_result(const std::string &dof_name) const
Get the result of a single DOF over time.
Definition Solver.cpp:316
void write_result_to_csv(const std::string &filename) const
Write the result to a csv file.
Definition Solver.cpp:383
std::string get_full_result() const
Get the full result as a csv encoded string.
Definition Solver.cpp:299
State of the system.
Definition State.h:46
DEBUG_MSG source file.
Simulation parameters.
Definition SimulationParameters.h:50