13#include <Kokkos_Core.hpp>
14#include <Kokkos_DualView.hpp>
16#if __has_include(<mkl_lapacke.h>)
17#include <mkl_lapacke.h>
22#include <KokkosBatched_Pbtrs.hpp>
23#include <KokkosBatched_Util.hpp>
25#include "splines_linear_problem.hpp"
27namespace ddc::detail {
30
31
32
33
34
35
36
37
38
39
40
41
42template <
class ExecSpace>
43class SplinesLinearProblemPDSBand :
public SplinesLinearProblem<ExecSpace>
46 using typename SplinesLinearProblem<ExecSpace>::MultiRHS;
47 using SplinesLinearProblem<ExecSpace>::size;
50 Kokkos::DualView<
double**, Kokkos::LayoutRight,
typename ExecSpace::memory_space>
55
56
57
58
59
60 explicit SplinesLinearProblemPDSBand(std::size_t
const mat_size, std::size_t
const kd)
61 : SplinesLinearProblem<ExecSpace>(mat_size)
62 , m_q(
"q", kd + 1, mat_size)
64 assert(m_q.extent(0) <= mat_size);
66 Kokkos::deep_copy(m_q.h_view, 0.);
69 double get_element(std::size_t i, std::size_t j)
const override
79 if (j - i < m_q.extent(0)) {
80 return m_q.h_view(j - i, i);
86 void set_element(std::size_t i, std::size_t j,
double const aij)
override
95 if (j - i < m_q.extent(0)) {
96 m_q.h_view(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_dpbtrf(
119 throw std::runtime_error(
120 "LAPACKE_dpbtrf failed with error code " + std::to_string(info));
129
130
131
132
133
134
135
136 void solve(MultiRHS
const b,
bool const)
const override
138 assert(b.extent(0) == size());
140 auto q_device = m_q.d_view;
141 Kokkos::RangePolicy<ExecSpace>
const policy(0, b.extent(1));
142 Kokkos::parallel_for(
145 KOKKOS_LAMBDA(
const int i) {
146 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
147 KokkosBatched::SerialPbtrs<
148 KokkosBatched::Uplo::Lower,
149 KokkosBatched::Algo::Pbtrs::Unblocked>::invoke(q_device, sub_b);
The top-level namespace of DDC.