BioCMAMC-ST
scalar_simulation.hpp
1#ifndef __SCALAR_SIMULATION_HPP__
2#define __SCALAR_SIMULATION_HPP__
3
4#include <Eigen/Core>
5#include <Eigen/Dense>
6#include <Eigen/Sparse>
7#include <Kokkos_Core.hpp>
8#include <cma_utils/cache_hydro_state.hpp>
9#include <common/common.hpp>
10#include <cstddef>
11#include <cstdint>
12#include <eigen_kokkos.hpp>
13#include <mc/traits.hpp>
14#include <simulation/alias.hpp>
15#include <simulation/mass_transfer.hpp>
16#include <span>
17#include <vector>
18
19namespace Simulation
20{
21
23 {
24 public:
25 ScalarSimulation(size_t n_compartments, size_t n_species, std::span<double> volume);
26
27 ScalarSimulation(ScalarSimulation&& other) noexcept = delete;
28 ScalarSimulation(const ScalarSimulation& other) noexcept = delete;
31 ~ScalarSimulation() = default;
32
33 bool deep_copy_concentration(const std::vector<double>& data);
34 void reduce_contribs(std::span<const double> data);
35
36 void performStepGL(double d_t,
37 const FlowMatrixType& m_transition,
38 const ColMajorMatrixtype& mtr,
40
41 void performStep(double d_t, const FlowMatrixType& m_transition);
42
43 // void performStep(double d_t, const FlowMatrixType& m_transition);
44
45 // Getters
46 [[nodiscard]] ColMajorMatrixtype& get_concentration();
47 [[nodiscard]] ColMajorKokkosScalarMatrix get_device_concentration() const;
48 [[nodiscard]] std::span<double const> getVolumeData() const;
49 [[nodiscard]] std::span<double> getContributionData() const;
50 [[nodiscard]] const DiagonalType& getVolume() const;
51 [[nodiscard]] auto getConcentrationArray() const;
52 [[nodiscard]] kernelContribution get_kernel_contribution() const;
53 [[nodiscard]] const ColMajorMatrixtype& get_mass_transfer() const;
54 [[nodiscard]] std::span<double> getConcentrationData();
55
56 [[nodiscard]] std::size_t n_row() const;
57 [[nodiscard]] std::size_t n_col() const;
58
59 // Setters
60
61 void set_mass();
62 void set_kernel_contribs_to_host() const;
63 void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val);
64 void set_sink(std::uint64_t i_compartment, double val);
65 void set_zero_contribs();
66 void setVolumes(std::span<const double> volumes, std::span<const double> inv_volumes);
67
68 private:
69 std::size_t n_r;
70 std::size_t n_c;
71 ColMajorMatrixtype total_mass;
72
73 DiagonalType volumes_inverse;
74 DiagonalType m_volumes;
75 DiagonalType sink;
78 };
79
81 {
82 // return alloc_concentrations.array();
83 return concentrations.eigen_data.array();
84 }
85
86 inline kernelContribution ScalarSimulation::get_kernel_contribution() const
87 {
88 return sources.compute;
89 }
90
92 {
93 sources.eigen_data.setZero();
94 Kokkos::deep_copy(sources.compute, 0);
95 }
96
101
102 inline void ScalarSimulation::set_feed(uint64_t i_r, uint64_t i_c, double val)
103 {
104 this->sources.eigen_data.coeffRef(EIGEN_INDEX(i_r), EIGEN_INDEX(i_c)) += val;
105 }
106
107 inline void ScalarSimulation::set_sink(uint64_t i_compartment, double val)
108 {
109 this->sink.diagonal().coeffRef(EIGEN_INDEX(i_compartment)) = val;
110 }
111
112 inline const DiagonalType& ScalarSimulation::getVolume() const
113 {
114 return m_volumes;
115 }
116
117 inline std::span<double> ScalarSimulation::getConcentrationData()
118 {
119 return this->concentrations.get_span();
120 }
121
122 inline std::span<double> ScalarSimulation::getContributionData() const
123 {
124 return {this->sources.host.data(), static_cast<size_t>(this->sources.host.size())};
125 }
126
127 inline std::span<double const> ScalarSimulation::getVolumeData() const
128 {
129 return {m_volumes.diagonal().data(), static_cast<size_t>(m_volumes.rows())};
130 }
131
132 inline void ScalarSimulation::setVolumes(std::span<const double> volumes,
133 std::span<const double> inv_volumes)
134 {
135 KOKKOS_ASSERT(volumes.size() == inv_volumes.size() && volumes.size() == n_col() &&
136 "scalar:setvolume")
137 // SIGFAULT ?
138 this->m_volumes.diagonal() =
139 Eigen::Map<const Eigen::VectorXd>(volumes.data(), static_cast<int>(volumes.size()));
140
141 this->volumes_inverse.diagonal() =
142 Eigen::Map<const Eigen::VectorXd>(inv_volumes.data(), static_cast<int>(inv_volumes.size()));
143 }
144
145 inline ScalarSimulation*
146 makeScalarSimulation(size_t n_compartments, size_t n_species, std::span<double> volumes)
147 {
148 return new ScalarSimulation(n_compartments, n_species, volumes); // NOLINT
149 }
150
151} // namespace Simulation
152
153#endif //__SCALAR_SIMULATION_HPP__
Definition scalar_simulation.hpp:23
ScalarSimulation(ScalarSimulation &&other) noexcept=delete
void performStepGL(double d_t, const FlowMatrixType &m_transition, const ColMajorMatrixtype &mtr, MassTransfer::Sign sign)
Definition implScalar.cpp:68
ScalarSimulation operator=(const ScalarSimulation &other)=delete
void set_sink(std::uint64_t i_compartment, double val)
Definition scalar_simulation.hpp:107
ColMajorKokkosScalarMatrix get_device_concentration() const
Definition implScalar.cpp:46
void setVolumes(std::span< const double > volumes, std::span< const double > inv_volumes)
Definition scalar_simulation.hpp:132
ScalarSimulation(const ScalarSimulation &other) noexcept=delete
std::size_t n_row() const
Definition implScalar.cpp:56
void set_feed(std::uint64_t i_r, std::uint64_t i_c, double val)
Definition scalar_simulation.hpp:102
std::size_t n_col() const
Definition implScalar.cpp:51
RowMajorEigenKokkos sources
Definition scalar_simulation.hpp:77
EigenKokkos concentrations
Definition scalar_simulation.hpp:76
void set_zero_contribs()
Definition scalar_simulation.hpp:91
ColMajorMatrixtype & get_concentration()
Definition implScalar.cpp:41
const DiagonalType & getVolume() const
Definition scalar_simulation.hpp:112
std::span< double const > getVolumeData() const
Definition scalar_simulation.hpp:127
DiagonalType volumes_inverse
Definition scalar_simulation.hpp:73
std::size_t n_r
Definition scalar_simulation.hpp:69
bool deep_copy_concentration(const std::vector< double > &data)
Definition implScalar.cpp:96
ScalarSimulation(size_t n_compartments, size_t n_species, std::span< double > volume)
Definition implScalar.cpp:14
void reduce_contribs(std::span< const double > data)
Definition implScalar.cpp:61
void performStep(double d_t, const FlowMatrixType &m_transition)
Definition implScalar.cpp:84
kernelContribution get_kernel_contribution() const
Definition scalar_simulation.hpp:86
ScalarSimulation operator=(ScalarSimulation &&other)=delete
auto getConcentrationArray() const
Definition scalar_simulation.hpp:80
std::span< double > getContributionData() const
Definition scalar_simulation.hpp:122
void set_mass()
Definition implScalar.cpp:108
const ColMajorMatrixtype & get_mass_transfer() const
ColMajorMatrixtype total_mass
Definition scalar_simulation.hpp:71
std::size_t n_c
Definition scalar_simulation.hpp:70
void set_kernel_contribs_to_host() const
Definition scalar_simulation.hpp:97
std::span< double > getConcentrationData()
Definition scalar_simulation.hpp:117
DiagonalType m_volumes
Definition scalar_simulation.hpp:74
DiagonalType sink
Definition scalar_simulation.hpp:75
Sign
Definition mass_transfer.hpp:37
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:146
ComputeView compute
Definition eigen_kokkos.hpp:98
EigenMatrix eigen_data
Definition eigen_kokkos.hpp:99
void update_compute_to_host() const
Definition eigen_kokkos.hpp:124
std::span< const double > get_span() const
Definition eigen_kokkos.hpp:109
HostView host
Definition eigen_kokkos.hpp:97