DDC 0.11.0
Loading...
Searching...
No Matches
splines_linear_problem_dense.cpp
1// Copyright (C) The DDC development team, see COPYRIGHT.md file
2//
3// SPDX-License-Identifier: MIT
4
5#include <cassert>
6#include <cstddef>
7#include <stdexcept>
8#include <string>
9
10#include <Kokkos_Core.hpp>
11
12#if __has_include(<mkl_lapacke.h>)
13# include <mkl_lapacke.h>
14#else
15# include <lapacke.h>
16#endif
17
18// The two following headers are necessary to workaround missing includes in KokkosBatched_Getrs.hpp.
19// They must be placed before `#include <KokkosBatched_Getrs.hpp>`.
20// This is fixed in Kokkos Kernels >=5.1.
21// clang-format off
22#include <KokkosBatched_Laswp.hpp>
23#include <KokkosBatched_Trsm_Decl.hpp>
24// clang-format on
25
26#include <KokkosBatched_Getrs.hpp>
27#include <KokkosBatched_Util.hpp>
28
31
32namespace ddc::detail {
33
34template <class ExecSpace>
35SplinesLinearProblemDense<ExecSpace>::SplinesLinearProblemDense(std::size_t const mat_size)
36 : SplinesLinearProblem<ExecSpace>(mat_size)
37 , m_a("a", mat_size, mat_size)
38 , m_ipiv("ipiv", mat_size)
39{
40 Kokkos::deep_copy(m_a.view_host(), 0.);
41}
42
43template <class ExecSpace>
44SplinesLinearProblemDense<ExecSpace>::~SplinesLinearProblemDense() = default;
45
46template <class ExecSpace>
47double SplinesLinearProblemDense<ExecSpace>::get_element(std::size_t const i, std::size_t const j)
48 const
49{
50 assert(i < size());
51 assert(j < size());
52 return m_a.view_host()(i, j);
53}
54
55template <class ExecSpace>
56void SplinesLinearProblemDense<ExecSpace>::set_element(
57 std::size_t const i,
58 std::size_t const j,
59 double const aij)
60{
61 assert(i < size());
62 assert(j < size());
63 m_a.view_host()(i, j) = aij;
64}
65
66template <class ExecSpace>
67void SplinesLinearProblemDense<ExecSpace>::setup_solver()
68{
69 int const info = LAPACKE_dgetrf(
70 LAPACK_ROW_MAJOR,
71 size(),
72 size(),
73 m_a.view_host().data(),
74 size(),
75 m_ipiv.view_host().data());
76 if (info != 0) {
77 throw std::runtime_error("LAPACKE_dgetrf failed with error code " + std::to_string(info));
78 }
79
80 // Convert 1-based index to 0-based index
81 for (std::size_t i = 0; i < size(); ++i) {
82 m_ipiv.view_host()(i) -= 1;
83 }
84
85 // Push on device
86 m_a.modify_host();
87 m_a.sync_device();
88 m_ipiv.modify_host();
89 m_ipiv.sync_device();
90}
91
92template <class ExecSpace>
93void SplinesLinearProblemDense<ExecSpace>::solve(MultiRHS const b, bool const transpose) const
94{
95 assert(b.extent(0) == size());
96
97 // For order 1 splines, size() can be 0 then we bypass the solver call.
98 if (size() == 0) {
99 return;
100 }
101
102 auto a_device = m_a.view_device();
103 auto ipiv_device = m_ipiv.view_device();
104
105 Kokkos::RangePolicy<ExecSpace> const policy(0, b.extent(1));
106
107 if (transpose) {
108 Kokkos::parallel_for(
109 "gerts",
110 policy,
111 KOKKOS_LAMBDA(int const i) {
112 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
113 KokkosBatched::SerialGetrs<
114 KokkosBatched::Trans::Transpose,
115 KokkosBatched::Algo::Getrs::Unblocked>::
116 invoke(a_device, ipiv_device, sub_b);
117 });
118 } else {
119 Kokkos::parallel_for(
120 "gerts",
121 policy,
122 KOKKOS_LAMBDA(int const i) {
123 auto sub_b = Kokkos::subview(b, Kokkos::ALL, i);
124 KokkosBatched::SerialGetrs<
125 KokkosBatched::Trans::NoTranspose,
126 KokkosBatched::Algo::Getrs::Unblocked>::
127 invoke(a_device, ipiv_device, sub_b);
128 });
129 }
130}
131
132#if defined(KOKKOS_ENABLE_SERIAL)
133template class SplinesLinearProblemDense<Kokkos::Serial>;
134#endif
135#if defined(KOKKOS_ENABLE_OPENMP)
136template class SplinesLinearProblemDense<Kokkos::OpenMP>;
137#endif
138#if defined(KOKKOS_ENABLE_CUDA)
139template class SplinesLinearProblemDense<Kokkos::Cuda>;
140#endif
141#if defined(KOKKOS_ENABLE_HIP)
142template class SplinesLinearProblemDense<Kokkos::HIP>;
143#endif
144#if defined(KOKKOS_ENABLE_SYCL)
145template class SplinesLinearProblemDense<Kokkos::SYCL>;
146#endif
147
148} // namespace ddc::detail
The top-level namespace of DDC.