16#include <Kokkos_Core.hpp>
17#include <Kokkos_DualView.hpp>
19#if __has_include(<mkl_lapacke.h>)
20# include <mkl_lapacke.h>
25#include <KokkosBatched_Pttrs.hpp>
26#include <KokkosBatched_Util.hpp>
28#include "splines_linear_problem.hpp"
30namespace ddc::detail {
33
34
35
36
37
38
39
40
41
42
43
44
45template <
class ExecSpace>
46class SplinesLinearProblemPDSTridiag :
public SplinesLinearProblem<ExecSpace>
49 using typename SplinesLinearProblem<ExecSpace>::MultiRHS;
50 using SplinesLinearProblem<ExecSpace>::size;
53 Kokkos::DualView<
double**, Kokkos::LayoutRight,
typename ExecSpace::memory_space>
58
59
60
61
62 explicit SplinesLinearProblemPDSTridiag(std::size_t
const mat_size)
63 : SplinesLinearProblem<ExecSpace>(mat_size)
64 , m_q(
"q", 2, mat_size)
66 Kokkos::deep_copy(m_q.view_host(), 0.);
69 double get_element(std::size_t i, std::size_t j)
const override
80 return m_q.view_host()(j - i, i);
86 void set_element(std::size_t i, std::size_t j,
double const aij)
override
96 m_q.view_host()(j - i, i) = aij;
98 assert(std::fabs(aij) < 1e-20);
103
104
105
106
107 void setup_solver()
override
109 int const info = LAPACKE_dpttrf(
111 m_q.view_host().data(),
112 m_q.view_host().data() + m_q.view_host().stride(0));
114 throw std::runtime_error(
115 "LAPACKE_dpttrf failed with error code " + std::to_string(info));
124
125
126
127
128
129
130
131 void solve(MultiRHS
const b,
bool const)
const override
133 assert(b.extent(0) == size());
134 auto q_device = m_q.view_device();
135 auto d = Kokkos::subview(q_device, 0, Kokkos::ALL);
137 subview(q_device, 1, Kokkos::pair<
int,
int>(0, q_device.extent_int(1) - 1));
138 Kokkos::RangePolicy<ExecSpace>
const policy(0, b.extent(1));
139 Kokkos::parallel_for(
142 KOKKOS_LAMBDA(
int const i) {
143 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
144 KokkosBatched::SerialPttrs<
145 KokkosBatched::Uplo::Lower,
146 KokkosBatched::Algo::Pttrs::Unblocked>::invoke(d, e, sub_b);
The top-level namespace of DDC.