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 <cma_utils/cache_hydro_state.hpp>
17#include <common/common.hpp>
18#include <cstddef>
19#include <cstdint>
20#include <eigen_kokkos.hpp>
21#include <mc/traits.hpp>
22#include <simulation/alias.hpp>
23#include <simulation/mass_transfer.hpp>
24#include <span>
25#include <vector>
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 void reduce_contribs(std::span<const double> data);
45
46 void performStepGL(double d_t,
47 const FlowMatrixType& m_transition,
48 const ColMajorMatrixtype& mtr,
50
51 void performStep(double d_t, const FlowMatrixType& m_transition);
52
53 void synchro_sources();
54
55 // void performStep(double d_t, const FlowMatrixType& m_transition);
56
57 // Getters
58 [[nodiscard]] ColMajorMatrixtype& get_concentration();
59 [[nodiscard]] ColMajorKokkosScalarMatrix get_device_concentration() const;
60 [[nodiscard]] std::span<double const> getVolumeData() const;
61 [[nodiscard]] std::span<double> getContributionData() const;
62
63 [[nodiscard]] std::span<double> getContributionData_mut();
64 [[nodiscard]] const DiagonalType& getVolume() const;
65 [[nodiscard]] auto getConcentrationArray() const;
66 [[nodiscard]] kernelContribution get_kernel_contribution() const;
67 [[nodiscard]] const ColMajorMatrixtype& get_mass_transfer() const;
68 [[nodiscard]] std::span<double> getConcentrationData();
69
70 [[nodiscard]] std::size_t n_row() const;
71 [[nodiscard]] std::size_t n_col() const;
72
73 // Setters
74
75 void set_mass();
76 void set_kernel_contribs_to_host() const;
77 void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val);
78 void set_sink(std::uint64_t i_compartment, double val);
79 void set_zero_contribs();
80 void setVolumes(std::span<const double> volumes,
81 std::span<const double> inv_volumes);
82
83 private:
84 std::size_t n_r;
85 std::size_t n_c;
86 ColMajorMatrixtype total_mass;
87 kernelContribution contribs;
88
89 DiagonalType volumes_inverse;
90 DiagonalType m_volumes;
91 DiagonalType sink;
92 EigenKokkos concentrations;
93 RowMajorEigenKokkos sources;
94 };
95
97 {
98 // return alloc_concentrations.array();
99 return concentrations.eigen_data.array();
100 }
101
102 inline kernelContribution ScalarSimulation::get_kernel_contribution() const
103 {
104 // return sources.compute;
105 return contribs;
106 }
107
109 {
110 sources.eigen_data.setZero();
111 Kokkos::deep_copy(contribs, 0);
112 Kokkos::deep_copy(sources.compute, 0);
113 this->sink.setZero();
114 }
115
117 {
118 sources.update_compute_to_host();
119 }
120
121 inline void ScalarSimulation::set_feed(uint64_t i_r, uint64_t i_c, double val)
122 {
123 this->sources.eigen_data.coeffRef(EIGEN_INDEX(i_r), EIGEN_INDEX(i_c)) +=
124 val;
125 }
126
127 inline void ScalarSimulation::set_sink(uint64_t i_compartment, double val)
128 {
129 this->sink.diagonal().coeffRef(EIGEN_INDEX(i_compartment)) += val;
130 }
131
132 inline const DiagonalType& ScalarSimulation::getVolume() const
133 {
134 return m_volumes;
135 }
136
137 inline std::span<double> ScalarSimulation::getConcentrationData()
138 {
139 return this->concentrations.get_span();
140 }
141
142 inline std::span<double> ScalarSimulation::getContributionData() const
143 {
144 return {this->sources.host.data(),
145 static_cast<size_t>(this->sources.host.size())};
146 }
147
149 {
150 return {this->sources.host.data(),
151 static_cast<size_t>(this->sources.host.size())};
152 }
153
154 inline std::span<double const> ScalarSimulation::getVolumeData() const
155 {
156 return {m_volumes.diagonal().data(), static_cast<size_t>(m_volumes.rows())};
157 }
158
159 inline void ScalarSimulation::setVolumes(std::span<const double> volumes,
160 std::span<const double> inv_volumes)
161 {
162 KOKKOS_ASSERT(volumes.size() == inv_volumes.size() &&
163 volumes.size() == n_col() && "scalar:setvolume")
164 // SIGFAULT ?
165 this->m_volumes.diagonal() = Eigen::Map<const Eigen::VectorXd>(
166 volumes.data(), static_cast<int>(volumes.size()));
167
168 this->volumes_inverse.diagonal() = Eigen::Map<const Eigen::VectorXd>(
169 inv_volumes.data(), static_cast<int>(inv_volumes.size()));
170 }
171
172 inline ScalarSimulation* makeScalarSimulation(size_t n_compartments,
173 size_t n_species,
174 std::span<double> volumes)
175 {
176 return new ScalarSimulation(n_compartments, n_species, volumes); // NOLINT
177 }
178
179} // namespace Simulation
180
181#endif //__SCALAR_SIMULATION_HPP__
Definition scalar_simulation.hpp:31
ScalarSimulation(ScalarSimulation &&other) noexcept=delete
void performStepGL(double d_t, const FlowMatrixType &m_transition, const ColMajorMatrixtype &mtr, MassTransfer::Sign sign)
Definition implScalar.cpp:84
ScalarSimulation operator=(const ScalarSimulation &other)=delete
void set_sink(std::uint64_t i_compartment, double val)
Definition scalar_simulation.hpp:127
ColMajorKokkosScalarMatrix get_device_concentration() const
Definition implScalar.cpp:61
void setVolumes(std::span< const double > volumes, std::span< const double > inv_volumes)
Definition scalar_simulation.hpp:159
std::span< double > getContributionData_mut()
Definition scalar_simulation.hpp:148
ScalarSimulation(const ScalarSimulation &other) noexcept=delete
std::size_t n_row() const
Definition implScalar.cpp:71
void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val)
Definition scalar_simulation.hpp:121
std::size_t n_col() const
Definition implScalar.cpp:66
RowMajorEigenKokkos sources
Definition scalar_simulation.hpp:93
EigenKokkos concentrations
Definition scalar_simulation.hpp:92
void set_zero_contribs()
Definition scalar_simulation.hpp:108
ColMajorMatrixtype & get_concentration()
Definition implScalar.cpp:55
kernelContribution contribs
Definition scalar_simulation.hpp:87
const DiagonalType & getVolume() const
Definition scalar_simulation.hpp:132
std::span< double const > getVolumeData() const
Definition scalar_simulation.hpp:154
DiagonalType volumes_inverse
Definition scalar_simulation.hpp:89
std::size_t n_r
Definition scalar_simulation.hpp:84
bool deep_copy_concentration(const std::vector< double > &data)
Definition implScalar.cpp:116
ScalarSimulation(size_t n_compartments, size_t n_species, std::span< double > volume)
Definition implScalar.cpp:21
void reduce_contribs(std::span< const double > data)
Definition implScalar.cpp:76
void performStep(double d_t, const FlowMatrixType &m_transition)
Definition implScalar.cpp:101
kernelContribution get_kernel_contribution() const
Definition scalar_simulation.hpp:102
ScalarSimulation operator=(ScalarSimulation &&other)=delete
auto getConcentrationArray() const
Definition scalar_simulation.hpp:96
void synchro_sources()
Definition implScalar.cpp:50
std::span< double > getContributionData() const
Definition scalar_simulation.hpp:142
void set_mass()
Definition implScalar.cpp:129
const ColMajorMatrixtype & get_mass_transfer() const
ColMajorMatrixtype total_mass
Definition scalar_simulation.hpp:86
std::size_t n_c
Definition scalar_simulation.hpp:85
void set_kernel_contribs_to_host() const
Definition scalar_simulation.hpp:116
std::span< double > getConcentrationData()
Definition scalar_simulation.hpp:137
DiagonalType m_volumes
Definition scalar_simulation.hpp:90
DiagonalType sink
Definition scalar_simulation.hpp:91
Sign
Definition mass_transfer.hpp:44
Namespace that contains classes and structures related to simulation handling.
Definition host_specific.hpp:13
ScalarSimulation * makeScalarSimulation(size_t n_compartments, size_t n_species, std::span< double > volumes)
Definition scalar_simulation.hpp:172