14#include <Kokkos_Core.hpp>
16#if __has_include
(<mkl_lapacke.h>)
17# include <mkl_lapacke.h>
22#include <KokkosBatched_Pttrs.hpp>
23#include <KokkosBatched_Util.hpp>
28namespace ddc::detail {
30template <
class ExecSpace>
31SplinesLinearProblemPDSTridiag<ExecSpace>::SplinesLinearProblemPDSTridiag(
32 std::size_t
const mat_size)
33 : SplinesLinearProblem<ExecSpace>(mat_size)
34 , m_q(
"q", 2, mat_size)
36 Kokkos::deep_copy(m_q.view_host(), 0.);
39template <
class ExecSpace>
40SplinesLinearProblemPDSTridiag<ExecSpace>::~SplinesLinearProblemPDSTridiag() =
default;
42template <
class ExecSpace>
43double SplinesLinearProblemPDSTridiag<ExecSpace>::get_element(std::size_t i, std::size_t j)
const
54 return m_q.view_host()(j - i, i);
60template <
class ExecSpace>
61void SplinesLinearProblemPDSTridiag<ExecSpace>::set_element(
74 m_q.view_host()(j - i, i) = aij;
76 assert(std::fabs(aij) < 1e-15);
80template <
class ExecSpace>
81void SplinesLinearProblemPDSTridiag<ExecSpace>::setup_solver()
83 int const info = LAPACKE_dpttrf(
85 m_q.view_host().data(),
86 m_q.view_host().data() + m_q.view_host().stride(0));
88 throw std::runtime_error(
"LAPACKE_dpttrf failed with error code " + std::to_string(info));
96template <
class ExecSpace>
97void SplinesLinearProblemPDSTridiag<ExecSpace>::solve(MultiRHS
const b,
bool const)
const
99 assert(b.extent(0) == size());
100 auto q_device = m_q.view_device();
101 auto d = Kokkos::subview(q_device, 0, Kokkos::ALL);
102 auto e = Kokkos::subview(q_device, 1, Kokkos::pair<
int,
int>(0, q_device.extent_int(1) - 1));
103 Kokkos::RangePolicy<ExecSpace>
const policy(0, b.extent(1));
104 Kokkos::parallel_for(
107 KOKKOS_LAMBDA(
int const i) {
108 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
109 KokkosBatched::SerialPttrs<
110 KokkosBatched::Uplo::Lower,
111 KokkosBatched::Algo::Pttrs::Unblocked>::invoke(d, e, sub_b);
115#if defined(KOKKOS_ENABLE_SERIAL
)
116template class SplinesLinearProblemPDSTridiag<Kokkos::Serial>;
118#if defined(KOKKOS_ENABLE_OPENMP)
119template class SplinesLinearProblemPDSTridiag<Kokkos::OpenMP>;
121#if defined(KOKKOS_ENABLE_CUDA)
122template class SplinesLinearProblemPDSTridiag<Kokkos::Cuda>;
124#if defined(KOKKOS_ENABLE_HIP)
125template class SplinesLinearProblemPDSTridiag<Kokkos::HIP>;
127#if defined(KOKKOS_ENABLE_SYCL)
128template class SplinesLinearProblemPDSTridiag<Kokkos::SYCL>;
The top-level namespace of DDC.