13#include <Kokkos_Core.hpp>
15#include "splines_linear_problem.hpp"
16#include "splines_linear_problem_2x2_blocks.hpp"
18namespace ddc::detail {
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37template <
class ExecSpace>
38class SplinesLinearProblem3x3Blocks :
public SplinesLinearProblem2x2Blocks<ExecSpace>
41 using typename SplinesLinearProblem2x2Blocks<ExecSpace>::MultiRHS;
42 using SplinesLinearProblem2x2Blocks<ExecSpace>::size;
43 using SplinesLinearProblem2x2Blocks<ExecSpace>::solve;
44 using SplinesLinearProblem2x2Blocks<ExecSpace>::m_top_left_block;
47 std::size_t m_top_size;
51
52
53
54
55
56
57 explicit SplinesLinearProblem3x3Blocks(
58 std::size_t
const mat_size,
59 std::size_t
const top_size,
60 std::unique_ptr<SplinesLinearProblem<ExecSpace>> center_block)
61 : SplinesLinearProblem2x2Blocks<ExecSpace>(mat_size, std::move(center_block))
62 , m_top_size(top_size)
68 void adjust_indices(std::size_t& i, std::size_t& j)
const
70 std::size_t
const nq = m_top_left_block->size();
74 }
else if (i < m_top_size + nq) {
80 }
else if (j < m_top_size + nq) {
86 double get_element(std::size_t i, std::size_t j)
const override
89 return SplinesLinearProblem2x2Blocks<ExecSpace>::get_element(i, j);
92 void set_element(std::size_t i, std::size_t j,
double const aij)
override
95 return SplinesLinearProblem2x2Blocks<ExecSpace>::set_element(i, j, aij);
100
101
102
103
104
105
106
107
108
109
110 void interchange_rows_from_3_to_2_blocks_rhs(MultiRHS
const b)
const
112 std::size_t
const nq = m_top_left_block->size();
114 MultiRHS
const b_top = Kokkos::
115 subview(b, std::pair<std::size_t, std::size_t> {0, m_top_size}, Kokkos::ALL);
116 MultiRHS
const b_bottom = Kokkos::
118 std::pair<std::size_t, std::size_t> {m_top_size + nq, size()},
121 MultiRHS
const b_top_dst = Kokkos::
123 std::pair<std::size_t, std::size_t> {m_top_size + nq, 2 * m_top_size + nq},
125 MultiRHS
const b_bottom_dst = Kokkos::
129 std::size_t> {2 * m_top_size + nq, m_top_size + size()},
132 if (b_bottom.extent(0) > b_top.extent(0)) {
134 MultiRHS
const buffer = Kokkos::create_mirror(ExecSpace(), b_bottom);
136 Kokkos::deep_copy(buffer, b_bottom);
137 Kokkos::deep_copy(b_bottom_dst, buffer);
139 Kokkos::deep_copy(b_bottom_dst, b_bottom);
141 Kokkos::deep_copy(b_top_dst, b_top);
145
146
147
148
149
150
151
152
153
154 void interchange_rows_from_2_to_3_blocks_rhs(MultiRHS
const b)
const
156 std::size_t
const nq = m_top_left_block->size();
158 MultiRHS
const b_top = Kokkos::
159 subview(b, std::pair<std::size_t, std::size_t> {0, m_top_size}, Kokkos::ALL);
160 MultiRHS
const b_bottom = Kokkos::
162 std::pair<std::size_t, std::size_t> {m_top_size + nq, size()},
165 MultiRHS
const b_top_src = Kokkos::
167 std::pair<std::size_t, std::size_t> {m_top_size + nq, 2 * m_top_size + nq},
169 MultiRHS
const b_bottom_src = Kokkos::
173 std::size_t> {2 * m_top_size + nq, m_top_size + size()},
176 Kokkos::deep_copy(b_top, b_top_src);
177 if (b_bottom.extent(0) > b_top.extent(0)) {
179 MultiRHS
const buffer = Kokkos::create_mirror(ExecSpace(), b_bottom);
181 Kokkos::deep_copy(buffer, b_bottom_src);
182 Kokkos::deep_copy(b_bottom, buffer);
184 Kokkos::deep_copy(b_bottom, b_bottom_src);
190
191
192
193
194
195
196
197
198
199 void solve(MultiRHS
const b,
bool const transpose)
const override
201 assert(b.extent(0) == size() + m_top_size);
203 interchange_rows_from_3_to_2_blocks_rhs(b);
204 SplinesLinearProblem2x2Blocks<ExecSpace>::
209 std::size_t> {m_top_size, m_top_size + size()},
212 interchange_rows_from_2_to_3_blocks_rhs(b);
216 std::size_t impl_required_number_of_rhs_rows()
const override
218 return size() + m_top_size;
The top-level namespace of DDC.