BioCMAMC-ST
scalar_simulation.hpp
1#ifndef __SCALAR_SIMULATION_HPP__
2#define __SCALAR_SIMULATION_HPP__
3
4#ifndef NDEBUG
5# pragma GCC diagnostic push
6# pragma GCC diagnostic ignored "-Wunused-but-set-variable"
7# pragma GCC diagnostic ignored "-Wnan-infinity-disabled"
8#endif
9#include <Eigen/Core>
10#include <Eigen/Dense>
11#include <Eigen/Sparse>
12#ifndef NDEBUG
13# pragma GCC diagnostic pop
14#endif
15#include <Kokkos_Core.hpp>
16#include <common/common.hpp>
17#include <cstddef>
18#include <cstdint>
19#include <eigen_kokkos.hpp>
20#include <mc/traits.hpp>
21#include <simulation/mass_transfer.hpp>
22#include <span>
23#include <vector>
24
25using FlowMatrixType = Eigen::SparseMatrix<double>;
26
27namespace Simulation
28{
29
31 {
32 public:
33 ScalarSimulation(size_t n_compartments,
34 size_t n_species,
35 std::span<double> volume);
36
37 ScalarSimulation(ScalarSimulation&& other) noexcept = delete;
38 ScalarSimulation(const ScalarSimulation& other) noexcept = delete;
41 ~ScalarSimulation() = default;
42
43 bool deep_copy_concentration(const std::vector<double>& data);
44
45 void reduce_contribs(std::span<const double> data);
46
48
49 void performStepGL(double d_t,
50 const ColMajorMatrixtype<double>& mtr,
52
53 void performStep(double d_t);
54
55 void synchro_sources();
56
57 // void performStep(double d_t, const FlowMatrixType& m_transition);
58
59 // Getters
60 [[nodiscard]] ColMajorMatrixtype<double>& get_concentration() noexcept;
61 [[nodiscard]] const DiagonalType& getVolume() const noexcept;
62 [[nodiscard]] std::size_t n_row() const noexcept;
63 [[nodiscard]] std::size_t n_col() const noexcept;
64
65 [[nodiscard]] ColMajorKokkosScalarMatrix<double>
67
68 [[nodiscard]] std::span<double const> volume_span() const;
69
70 [[nodiscard]] std::span<double> contribution_span() const;
71
72 [[nodiscard]] std::span<double> contribution_span_mut();
73
74 [[nodiscard]] auto getConcentrationArray() const;
75 [[nodiscard]] kernelContribution get_kernel_contribution() const;
76
77 [[nodiscard]] std::span<double> getConcentrationData();
78
79 [[nodiscard]] std::span<const double> getConcentrationData() const;
80
81 // Setters
82
83 void set_mass();
84 void set_kernel_contribs_to_host() const;
85 void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val);
86 void set_sink(std::uint64_t i_compartment, double val);
87 void set_zero_contribs();
88 void setVolumes(std::span<const double> volumes,
89 std::span<const double> inv_volumes);
90
91 private:
92 std::size_t n_r;
93 std::size_t n_c;
94 ColMajorMatrixtype<double> total_mass;
95 kernelContribution contribs;
96
97 FlowMatrixType m_transition;
98
99 DiagonalType volumes_inverse;
100 DiagonalType m_volumes;
101 DiagonalType sink;
102 EigenKokkos<double> concentrations;
103 RowMajorEigenKokkos<double> sources;
104 };
105
106 inline auto
108 {
109 // return alloc_concentrations.array();
110 return concentrations.eigen_data.array();
111 }
112
113 inline kernelContribution
115 {
116 // return sources.compute;
117 return contribs;
118 }
119
120 inline void
122 {
123 sources.eigen_data.setZero();
124 Kokkos::deep_copy(contribs, 0);
125 Kokkos::deep_copy(sources.compute, 0);
126 this->sink.setZero();
127 }
128
129 inline void
131 {
132 sources.update_compute_to_host();
133 }
134
135 inline void
136 ScalarSimulation::set_feed(uint64_t i_r, uint64_t i_c, double val)
137 {
138 this->sources.eigen_data.coeffRef(EIGEN_INDEX(i_r), EIGEN_INDEX(i_c))
139 += val;
140 }
141
142 inline void
143 ScalarSimulation::set_sink(uint64_t i_compartment, double val)
144 {
145 this->sink.diagonal().coeffRef(EIGEN_INDEX(i_compartment)) += val;
146 }
147
148 inline const DiagonalType&
150 {
151 return m_volumes;
152 }
153
154 inline std::span<double>
156 {
157 return this->concentrations.get_span();
158 }
159
160 inline std::span<const double>
162 {
163 return this->concentrations.get_span();
164 }
165
166 inline std::span<double>
168 {
169 return { this->sources.host.data(),
170 static_cast<size_t>(this->sources.host.size()) };
171 }
172
173 inline std::span<double>
175 {
176 return { this->sources.host.data(),
177 static_cast<size_t>(this->sources.host.size()) };
178 }
179
180 inline std::span<double const>
182 {
183 return { m_volumes.diagonal().data(),
184 static_cast<size_t>(m_volumes.rows()) };
185 }
186
187 inline void
188 ScalarSimulation::setVolumes(std::span<const double> volumes,
189 std::span<const double> inv_volumes)
190 {
191 KOKKOS_ASSERT(volumes.size() == inv_volumes.size()
192 && volumes.size() == n_col() && "scalar:setvolume")
193 // SIGFAULT ?
194 this->m_volumes.diagonal() = Eigen::Map<const Eigen::VectorXd>(
195 volumes.data(), static_cast<int>(volumes.size()));
196
197 this->volumes_inverse.diagonal() = Eigen::Map<const Eigen::VectorXd>(
198 inv_volumes.data(), static_cast<int>(inv_volumes.size()));
199 }
200
201 inline ScalarSimulation*
202 makeScalarSimulation(size_t n_compartments,
203 size_t n_species,
204 std::span<double> volumes)
205 {
206 return new ScalarSimulation(n_compartments, n_species, volumes); // NOLINT
207 }
208
209} // namespace Simulation
210
211#endif //__SCALAR_SIMULATION_HPP__
Definition scalar_simulation.hpp:31
ScalarSimulation(ScalarSimulation &&other) noexcept=delete
ScalarSimulation operator=(const ScalarSimulation &other)=delete
void performStepGL(double d_t, const ColMajorMatrixtype< double > &mtr, MassTransfer::Sign sign)
Definition implScalar.cpp:178
void performStep(double d_t)
Definition implScalar.cpp:195
void set_sink(std::uint64_t i_compartment, double val)
Definition scalar_simulation.hpp:143
ColMajorMatrixtype< double > & get_concentration() noexcept
Definition implScalar.cpp:139
void setVolumes(std::span< const double > volumes, std::span< const double > inv_volumes)
Definition scalar_simulation.hpp:188
ScalarSimulation(const ScalarSimulation &other) noexcept=delete
void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val)
Definition scalar_simulation.hpp:136
std::size_t n_col() const noexcept
Definition implScalar.cpp:145
EigenKokkos< double > concentrations
Definition scalar_simulation.hpp:102
std::span< double > contribution_span() const
Definition scalar_simulation.hpp:167
void set_zero_contribs()
Definition scalar_simulation.hpp:121
kernelContribution contribs
Definition scalar_simulation.hpp:95
std::span< double > contribution_span_mut()
Definition scalar_simulation.hpp:174
void set_transition(CmaUtils::StateCooMatrixType &&transition)
Definition implScalar.cpp:229
DiagonalType volumes_inverse
Definition scalar_simulation.hpp:99
std::size_t n_r
Definition scalar_simulation.hpp:92
bool deep_copy_concentration(const std::vector< double > &data)
Definition implScalar.cpp:209
ScalarSimulation(size_t n_compartments, size_t n_species, std::span< double > volume)
Definition implScalar.cpp:103
void reduce_contribs(std::span< const double > data)
Definition implScalar.cpp:169
std::size_t n_row() const noexcept
Definition implScalar.cpp:151
FlowMatrixType m_transition
Definition scalar_simulation.hpp:97
kernelContribution get_kernel_contribution() const
Definition scalar_simulation.hpp:114
const DiagonalType & getVolume() const noexcept
Definition scalar_simulation.hpp:149
std::span< double const > volume_span() const
Definition scalar_simulation.hpp:181
ColMajorKokkosScalarMatrix< double > get_device_concentration() const
Definition implScalar.cpp:163
ScalarSimulation operator=(ScalarSimulation &&other)=delete
auto getConcentrationArray() const
Definition scalar_simulation.hpp:107
void synchro_sources()
Definition implScalar.cpp:157
void set_mass()
Definition implScalar.cpp:223
RowMajorEigenKokkos< double > sources
Definition scalar_simulation.hpp:103
ColMajorMatrixtype< double > total_mass
Definition scalar_simulation.hpp:94
std::size_t n_c
Definition scalar_simulation.hpp:93
void set_kernel_contribs_to_host() const
Definition scalar_simulation.hpp:130
std::span< double > getConcentrationData()
Definition scalar_simulation.hpp:155
DiagonalType m_volumes
Definition scalar_simulation.hpp:100
DiagonalType sink
Definition scalar_simulation.hpp:101
::rust::Box<::CooMatrixWrap > StateCooMatrixType
Opaque type for Naive COO matrix.
Definition alias.hpp:26
Sign
Definition mass_transfer.hpp:44
Namespace that contains classes and structures related to simulation handling.
Definition host_specific.hpp:12
ScalarSimulation * makeScalarSimulation(size_t n_compartments, size_t n_species, std::span< double > volumes)
Definition scalar_simulation.hpp:202