BioCMAMC-ST
cache_hydro_state.hpp
1#ifndef __CACHE_HYDRO_STATE__
2#define __CACHE_HYDRO_STATE__
3
4#include <Eigen/Core>
5#include <Eigen/Dense>
6#include <Eigen/Sparse>
7#include <common/common.hpp>
8#include <vector>
9
10template <typename ExecSpace>
11using DiagonalView = Kokkos::
12 View<double*, Kokkos::LayoutLeft, ExecSpace, Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
13
14template <typename Space>
15using CumulativeProbabilityView =
16 Kokkos::View<double**, Kokkos::LayoutRight, Space, Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
17
18using FlowMatrixType = Eigen::SparseMatrix<double>;
19
26namespace CmaUtils
27{
28 class ProxyPreCalculatedHydroState;
29
40 {
41 public:
43
48
53
58 explicit PreCalculatedHydroState(const FlowMatrixType& _tm);
59
80 [[nodiscard]] const FlowMatrixType& get_transition() const;
81
86 [[nodiscard]] DiagonalView<ComputeSpace> get_kernel_diagonal() const;
87
88 FlowMatrixType transition_matrix;
89
90 Eigen::Matrix<double, -1, -1, Eigen::RowMajor>
92
93 std::vector<double> inverse_volume;
94 std::vector<double> volume;
95 DiagonalView<ComputeSpace> diagonal_compute;
96
97 // CumulativeProbabilityView<ComputeSpace>
98 // compute_cumulative_probability; ///< View for cumulative probability computation.
99
100 private:
101 };
102
103 inline DiagonalView<ComputeSpace> PreCalculatedHydroState::get_kernel_diagonal() const
104 {
105 return diagonal_compute;
106 }
107
108} // namespace CmaUtils
109
110#endif //__CACHE_HYDRO_STATE__
Structure to store hydrodynamic properties from a compartment mesh.
Definition cache_hydro_state.hpp:40
PreCalculatedHydroState(const PreCalculatedHydroState &other)=delete
Deleted copy assigment.
FlowMatrixType transition_matrix
Matrix representing flow transitions.
Definition cache_hydro_state.hpp:88
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:91
std::vector< double > inverse_volume
Inverse of compartment volumes.
Definition cache_hydro_state.hpp:93
PreCalculatedHydroState()
Default constructor.
Definition cache_hydro_state.cpp:15
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:94
friend class ProxyCalculatedHydroState
Definition cache_hydro_state.hpp:42
DiagonalView< ComputeSpace > diagonal_compute
Diagonal view for compute operations.
Definition cache_hydro_state.hpp:95
DiagonalView< ComputeSpace > get_kernel_diagonal() const
Get the kernel diagonal view.
Definition cache_hydro_state.hpp:103
const FlowMatrixType & get_transition() const
Get the transition matrix.
Definition cache_hydro_state.cpp:20
Namespace to handle algorithms and structures related to reading compartment mesh.
Definition host_specific.hpp:17