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_Pbtrs.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 SplinesLinearProblemPDSBand :
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
63 explicit SplinesLinearProblemPDSBand(std::size_t
const mat_size, std::size_t
const kd)
64 : SplinesLinearProblem<ExecSpace>(mat_size)
65 , m_q(
"q", kd + 1, mat_size)
67 assert(m_q.extent(0) <= mat_size);
69 Kokkos::deep_copy(m_q.view_host(), 0.);
72 double get_element(std::size_t i, std::size_t j)
const override
82 if (j - i < m_q.extent(0)) {
83 return m_q.view_host()(j - i, i);
89 void set_element(std::size_t i, std::size_t j,
double const aij)
override
98 if (j - i < m_q.extent(0)) {
99 m_q.view_host()(j - i, i) = aij;
101 assert(std::fabs(aij) < 1e-20);
106
107
108
109
110 void setup_solver()
override
112 int const info = LAPACKE_dpbtrf(
117 m_q.view_host().data(),
118 m_q.view_host().stride(
122 throw std::runtime_error(
123 "LAPACKE_dpbtrf failed with error code " + std::to_string(info));
132
133
134
135
136
137
138
139 void solve(MultiRHS
const b,
bool const)
const override
141 assert(b.extent(0) == size());
143 auto q_device = m_q.view_device();
144 Kokkos::RangePolicy<ExecSpace>
const policy(0, b.extent(1));
145 Kokkos::parallel_for(
148 KOKKOS_LAMBDA(
int const i) {
149 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
150 KokkosBatched::SerialPbtrs<
151 KokkosBatched::Uplo::Lower,
152 KokkosBatched::Algo::Pbtrs::Unblocked>::invoke(q_device, sub_b);
The top-level namespace of DDC.