BioCMAMC-ST
cache_hydro_state.hpp
1#ifndef __CACHE_HYDRO_STATE__
2#define __CACHE_HYDRO_STATE__
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
16#include <common/common.hpp>
17#include <vector>
18
19template <typename ExecSpace>
20using DiagonalView = Kokkos::View<double*,
21 Kokkos::LayoutLeft,
22 ExecSpace,
23 Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
24
25template <typename Space>
26using CumulativeProbabilityView =
27 Kokkos::View<double**,
28 Kokkos::LayoutRight,
29 Space,
30 Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
31
32using FlowMatrixType = Eigen::SparseMatrix<double>;
33
41namespace CmaUtils
42{
44
56 {
57 public:
59
64
69
74 explicit PreCalculatedHydroState(const FlowMatrixType& _tm);
75
84 operator=(const PreCalculatedHydroState& rhs) = delete;
97 [[nodiscard]] const FlowMatrixType& get_transition() const;
98
103 [[nodiscard]] DiagonalView<ComputeSpace> get_kernel_diagonal() const;
104
105 FlowMatrixType transition_matrix;
106
107 Eigen::Matrix<double, -1, -1, Eigen::RowMajor>
109
110 std::vector<double> inverse_volume;
111 std::vector<double> volume;
112 DiagonalView<ComputeSpace>
114
115 // CumulativeProbabilityView<ComputeSpace>
116 // compute_cumulative_probability; ///< View for cumulative probability
117 // computation.
118
119 private:
120 };
121
122 inline DiagonalView<ComputeSpace>
127
128} // namespace CmaUtils
129
130#endif //__CACHE_HYDRO_STATE__
PreCalculatedHydroState(const PreCalculatedHydroState &other)=delete
Deleted copy assigment.
FlowMatrixType transition_matrix
Matrix representing flow transitions.
Definition cache_hydro_state.hpp:105
PreCalculatedHydroState & operator=(const PreCalculatedHydroState &rhs)=delete
Deleted copy constructor.
Eigen::Matrix< double, -1, -1, Eigen::RowMajor > cumulative_probability
Cumulative probability matrix.
Definition cache_hydro_state.hpp:108
std::vector< double > inverse_volume
Inverse of compartment volumes.
Definition cache_hydro_state.hpp:110
PreCalculatedHydroState()
Default constructor.
Definition cache_hydro_state.cpp:13
PreCalculatedHydroState(PreCalculatedHydroState &&other)=default
Default move assigment.
PreCalculatedHydroState & operator=(PreCalculatedHydroState &&rhs)=default
Default move constructor.
~PreCalculatedHydroState()=default
Default destructor.
std::vector< double > volume
Volumes of compartments.
Definition cache_hydro_state.hpp:111
friend class ProxyCalculatedHydroState
Definition cache_hydro_state.hpp:58
DiagonalView< ComputeSpace > diagonal_compute
Diagonal view for compute operations.
Definition cache_hydro_state.hpp:113
DiagonalView< ComputeSpace > get_kernel_diagonal() const
Get the kernel diagonal view.
Definition cache_hydro_state.hpp:123
const FlowMatrixType & get_transition() const
Get the transition matrix.
Definition cache_hydro_state.cpp:20
Proxy class for PreCalculatedHydroState used to fill this struct with data from an external library.
Definition proxy_cache.hpp:28
Namespace to handle algorithms and structures related to reading compartment mesh.
Definition host_specific.hpp:18