10#ifndef IFPACK2_BLOCKCOMPUTERES_IMPL_HPP
11#define IFPACK2_BLOCKCOMPUTERES_IMPL_HPP
13#include "Ifpack2_BlockHelper.hpp"
17namespace BlockHelperDetails {
22template <
typename MatrixType>
25 using local_ordinal_type_1d_view =
typename impl_type::local_ordinal_type_1d_view;
27 using i64_3d_view =
typename impl_type::i64_3d_view;
28 using impl_scalar_type_1d_view_tpetra = Unmanaged<typename impl_type::impl_scalar_type_1d_view_tpetra>;
30 size_type_1d_view rowptr, rowptr_remote;
37 local_ordinal_type_1d_view A_colindsub, A_colindsub_remote;
39 i64_3d_view A_x_offsets;
41 i64_3d_view A_x_offsets_remote;
44 bool is_tpetra_block_crs;
47 impl_scalar_type_1d_view_tpetra tpetra_values;
50 AmD(
const AmD &b) =
default;
53template <
typename MatrixType>
55 using local_ordinal_type =
typename BlockHelperDetails::ImplType<MatrixType>::local_ordinal_type;
56 using local_ordinal_type_1d_view =
typename BlockHelperDetails::ImplType<MatrixType>::local_ordinal_type_1d_view;
57 using local_ordinal_type_2d_view =
typename BlockHelperDetails::ImplType<MatrixType>::local_ordinal_type_2d_view;
59 PartInterface() =
default;
60 PartInterface(
const PartInterface &b) =
default;
80 local_ordinal_type_1d_view lclrow;
82 local_ordinal_type_1d_view partptr;
83 local_ordinal_type_2d_view partptr_sub;
84 local_ordinal_type_1d_view partptr_schur;
87 local_ordinal_type_1d_view packptr;
88 local_ordinal_type_1d_view packptr_sub;
89 local_ordinal_type_1d_view packindices_sub;
90 local_ordinal_type_2d_view packindices_schur;
93 local_ordinal_type_1d_view part2rowidx0;
94 local_ordinal_type_1d_view part2rowidx0_sub;
98 local_ordinal_type_1d_view part2packrowidx0;
99 local_ordinal_type_2d_view part2packrowidx0_sub;
100 local_ordinal_type part2packrowidx0_back;
102 local_ordinal_type_1d_view rowidx2part;
103 local_ordinal_type_1d_view rowidx2part_sub;
111 local_ordinal_type max_partsz;
112 local_ordinal_type max_subpartsz;
113 local_ordinal_type n_subparts_per_part;
114 local_ordinal_type nparts;
124template <
typename MatrixType>
125void precompute_A_x_offsets(
126 AmD<MatrixType> &amd,
127 const PartInterface<MatrixType> &interf,
128 const Teuchos::RCP<
const typename ImplType<MatrixType>::tpetra_crs_graph_type> &g,
129 const typename ImplType<MatrixType>::local_ordinal_type_1d_view &dm2cm,
131 bool ownedRemoteSeparate) {
132 using impl_type = ImplType<MatrixType>;
133 using i64_3d_view =
typename impl_type::i64_3d_view;
134 using size_type =
typename impl_type::size_type;
135 using local_ordinal_type =
typename impl_type::local_ordinal_type;
136 using execution_space =
typename impl_type::execution_space;
137 auto local_graph = g->getLocalGraphDevice();
138 const auto A_block_rowptr = local_graph.row_map;
139 const auto A_colind = local_graph.entries;
140 local_ordinal_type numLocalRows = interf.rowidx2part.extent(0);
141 int blocksize_square = blocksize * blocksize;
143 auto lclrow = interf.lclrow;
144 auto A_colindsub = amd.A_colindsub;
145 auto A_colindsub_remote = amd.A_colindsub_remote;
146 auto rowptr = amd.rowptr;
147 auto rowptr_remote = amd.rowptr_remote;
148 bool is_dm2cm_active = dm2cm.extent(0);
149 if (ownedRemoteSeparate) {
151 local_ordinal_type maxOwnedEntriesPerRow = 0;
152 local_ordinal_type maxNonownedEntriesPerRow = 0;
153 Kokkos::parallel_reduce(
154 Kokkos::RangePolicy<execution_space>(0, numLocalRows),
155 KOKKOS_LAMBDA(local_ordinal_type i, local_ordinal_type & lmaxOwned, local_ordinal_type & lmaxNonowned) {
156 const local_ordinal_type lr = lclrow(i);
157 local_ordinal_type rowNumOwned = rowptr(lr + 1) - rowptr(lr);
158 if (rowNumOwned > lmaxOwned)
159 lmaxOwned = rowNumOwned;
161 if (rowptr_remote.extent(0)) {
162 local_ordinal_type rowNumNonowned = rowptr_remote(lr + 1) - rowptr_remote(lr);
163 if (rowNumNonowned > lmaxNonowned)
164 lmaxNonowned = rowNumNonowned;
169 Kokkos::Max<local_ordinal_type>(maxOwnedEntriesPerRow), Kokkos::Max<local_ordinal_type>(maxNonownedEntriesPerRow));
173 amd.A_x_offsets = i64_3d_view(
"amd.A_x_offsets", numLocalRows, 2, maxOwnedEntriesPerRow);
174 amd.A_x_offsets_remote = i64_3d_view(
"amd.A_x_offsets_remote", numLocalRows, 2, maxNonownedEntriesPerRow);
175 auto A_x_offsets = amd.A_x_offsets;
176 auto A_x_offsets_remote = amd.A_x_offsets_remote;
178 Kokkos::parallel_for(
179 Kokkos::RangePolicy<execution_space>(0, numLocalRows),
180 KOKKOS_LAMBDA(local_ordinal_type i) {
181 const local_ordinal_type lr = lclrow(i);
182 const size_type A_k0 = A_block_rowptr(lr);
184 size_type rowBegin = rowptr(lr);
185 local_ordinal_type rowNumOwned = rowptr(lr + 1) - rowBegin;
186 for (local_ordinal_type entry = 0; entry < maxOwnedEntriesPerRow; entry++) {
187 if (entry < rowNumOwned) {
188 const size_type j = A_k0 + A_colindsub(rowBegin + entry);
189 const local_ordinal_type A_colind_at_j = A_colind(j);
190 const local_ordinal_type loc = is_dm2cm_active ? dm2cm(A_colind_at_j) : A_colind_at_j;
191 A_x_offsets(i, 0, entry) = int64_t(j) * blocksize_square;
192 A_x_offsets(i, 1, entry) = int64_t(loc) * blocksize;
194 A_x_offsets(i, 0, entry) = KokkosKernels::ArithTraits<int64_t>::min();
195 A_x_offsets(i, 1, entry) = KokkosKernels::ArithTraits<int64_t>::min();
199 if (rowptr_remote.extent(0)) {
200 rowBegin = rowptr_remote(lr);
201 local_ordinal_type rowNumNonowned = rowptr_remote(lr + 1) - rowBegin;
202 for (local_ordinal_type entry = 0; entry < maxNonownedEntriesPerRow; entry++) {
203 if (entry < rowNumNonowned) {
204 const size_type j = A_k0 + A_colindsub_remote(rowBegin + entry);
205 const local_ordinal_type A_colind_at_j = A_colind(j);
206 const local_ordinal_type loc = A_colind_at_j - numLocalRows;
207 A_x_offsets_remote(i, 0, entry) = int64_t(j) * blocksize_square;
208 A_x_offsets_remote(i, 1, entry) = int64_t(loc) * blocksize;
210 A_x_offsets_remote(i, 0, entry) = KokkosKernels::ArithTraits<int64_t>::min();
211 A_x_offsets_remote(i, 1, entry) = KokkosKernels::ArithTraits<int64_t>::min();
218 local_ordinal_type maxEntriesPerRow = 0;
219 Kokkos::parallel_reduce(
220 Kokkos::RangePolicy<execution_space>(0, numLocalRows),
221 KOKKOS_LAMBDA(local_ordinal_type i, local_ordinal_type & lmax) {
222 const local_ordinal_type lr = lclrow(i);
223 local_ordinal_type rowNum = rowptr(lr + 1) - rowptr(lr);
227 Kokkos::Max<local_ordinal_type>(maxEntriesPerRow));
228 amd.A_x_offsets = i64_3d_view(
"amd.A_x_offsets", numLocalRows, 2, maxEntriesPerRow);
229 auto A_x_offsets = amd.A_x_offsets;
232 Kokkos::parallel_for(
233 Kokkos::RangePolicy<execution_space>(0, numLocalRows),
234 KOKKOS_LAMBDA(local_ordinal_type i) {
235 const local_ordinal_type lr = lclrow(i);
236 const size_type A_k0 = A_block_rowptr(lr);
238 size_type rowBegin = rowptr(lr);
239 local_ordinal_type rowOwned = rowptr(lr + 1) - rowBegin;
240 for (local_ordinal_type entry = 0; entry < maxEntriesPerRow; entry++) {
241 if (entry < rowOwned) {
242 const size_type j = A_k0 + A_colindsub(rowBegin + entry);
243 A_x_offsets(i, 0, entry) = j * blocksize_square;
244 const local_ordinal_type A_colind_at_j = A_colind(j);
245 if (A_colind_at_j < numLocalRows) {
246 const local_ordinal_type loc = is_dm2cm_active ? dm2cm[A_colind_at_j] : A_colind_at_j;
247 A_x_offsets(i, 1, entry) = int64_t(loc) * blocksize;
249 A_x_offsets(i, 1, entry) = int64_t(A_colind_at_j) * blocksize;
252 A_x_offsets(i, 0, entry) = KokkosKernels::ArithTraits<int64_t>::min();
253 A_x_offsets(i, 1, entry) = KokkosKernels::ArithTraits<int64_t>::min();
263static inline int ComputeResidualVectorRecommendedCudaVectorSize(
const int blksize,
264 const int team_size) {
265 int total_team_size(0);
267 total_team_size = 32;
268 else if (blksize <= 9)
269 total_team_size = 32;
270 else if (blksize <= 12)
271 total_team_size = 96;
272 else if (blksize <= 16)
273 total_team_size = 128;
274 else if (blksize <= 20)
275 total_team_size = 160;
277 total_team_size = 160;
278 return total_team_size / team_size;
281static inline int ComputeResidualVectorRecommendedHIPVectorSize(
const int blksize,
282 const int team_size) {
283 int total_team_size(0);
285 total_team_size = 32;
286 else if (blksize <= 9)
287 total_team_size = 32;
288 else if (blksize <= 12)
289 total_team_size = 96;
290 else if (blksize <= 16)
291 total_team_size = 128;
292 else if (blksize <= 20)
293 total_team_size = 160;
295 total_team_size = 160;
296 return total_team_size / team_size;
299static inline int ComputeResidualVectorRecommendedSYCLVectorSize(
const int blksize,
300 const int team_size) {
301 int total_team_size(0);
303 total_team_size = 32;
304 else if (blksize <= 9)
305 total_team_size = 32;
306 else if (blksize <= 12)
307 total_team_size = 96;
308 else if (blksize <= 16)
309 total_team_size = 128;
310 else if (blksize <= 20)
311 total_team_size = 160;
313 total_team_size = 160;
314 return total_team_size / team_size;
318static inline int ComputeResidualVectorRecommendedVectorSize(
const int blksize,
319 const int team_size) {
320 if (is_cuda<T>::value)
321 return ComputeResidualVectorRecommendedCudaVectorSize(blksize, team_size);
322 if (is_hip<T>::value)
323 return ComputeResidualVectorRecommendedHIPVectorSize(blksize, team_size);
324 if (is_sycl<T>::value)
325 return ComputeResidualVectorRecommendedSYCLVectorSize(blksize, team_size);
329template <
typename MatrixType>
330struct ComputeResidualVector {
332 using impl_type = BlockHelperDetails::ImplType<MatrixType>;
334 using execution_space =
typename impl_type::execution_space;
335 using memory_space =
typename impl_type::memory_space;
337 using local_ordinal_type =
typename impl_type::local_ordinal_type;
340 using magnitude_type =
typename impl_type::magnitude_type;
341 using btdm_scalar_type =
typename impl_type::btdm_scalar_type;
342 using btdm_magnitude_type =
typename impl_type::btdm_magnitude_type;
344 using local_ordinal_type_1d_view =
typename impl_type::local_ordinal_type_1d_view;
346 using tpetra_block_access_view_type =
typename impl_type::tpetra_block_access_view_type;
347 using impl_scalar_type_1d_view =
typename impl_type::impl_scalar_type_1d_view;
348 using impl_scalar_type_2d_view_tpetra =
typename impl_type::impl_scalar_type_2d_view_tpetra;
349 using vector_type_3d_view =
typename impl_type::vector_type_3d_view;
350 using btdm_scalar_type_4d_view =
typename impl_type::btdm_scalar_type_4d_view;
351 using i64_3d_view =
typename impl_type::i64_3d_view;
352 static constexpr int vector_length = impl_type::vector_length;
355 using member_type =
typename Kokkos::TeamPolicy<execution_space>::member_type;
358 enum :
int { max_blocksize = 32 };
361 ConstUnmanaged<impl_scalar_type_2d_view_tpetra> b;
362 ConstUnmanaged<impl_scalar_type_2d_view_tpetra> x;
363 ConstUnmanaged<impl_scalar_type_2d_view_tpetra> x_remote;
364 Unmanaged<impl_scalar_type_2d_view_tpetra> y;
365 Unmanaged<vector_type_3d_view> y_packed;
366 Unmanaged<btdm_scalar_type_4d_view> y_packed_scalar;
369 const ConstUnmanaged<size_type_1d_view> rowptr, rowptr_remote;
370 const ConstUnmanaged<local_ordinal_type_1d_view> colindsub, colindsub_remote;
371 const ConstUnmanaged<impl_scalar_type_1d_view> tpetra_values;
375 const ConstUnmanaged<Kokkos::View<size_t *, node_device_type>> A_block_rowptr;
376 const ConstUnmanaged<Kokkos::View<size_t *, node_device_type>> A_point_rowptr;
377 const ConstUnmanaged<Kokkos::View<local_ordinal_type *, node_device_type>> A_colind;
380 const local_ordinal_type blocksize_requested;
383 const ConstUnmanaged<local_ordinal_type_1d_view> part2packrowidx0;
384 const ConstUnmanaged<local_ordinal_type_1d_view> part2rowidx0;
385 const ConstUnmanaged<local_ordinal_type_1d_view> rowidx2part;
386 const ConstUnmanaged<local_ordinal_type_1d_view> partptr;
387 const ConstUnmanaged<local_ordinal_type_1d_view> lclrow;
388 const ConstUnmanaged<local_ordinal_type_1d_view> dm2cm;
391 const ConstUnmanaged<i64_3d_view> A_x_offsets;
392 const ConstUnmanaged<i64_3d_view> A_x_offsets_remote;
394 const bool is_dm2cm_active;
395 const bool hasBlockCrsMatrix;
398 template <
typename LocalCrsGraphType>
399 ComputeResidualVector(
const AmD<MatrixType> &amd,
400 const LocalCrsGraphType &block_graph,
401 const LocalCrsGraphType &point_graph,
402 const local_ordinal_type &blocksize_requested_,
403 const PartInterface<MatrixType> &interf,
404 const local_ordinal_type_1d_view &dm2cm_,
405 bool hasBlockCrsMatrix_)
407 , rowptr_remote(amd.rowptr_remote)
408 , colindsub(amd.A_colindsub)
409 , colindsub_remote(amd.A_colindsub_remote)
410 , tpetra_values(amd.tpetra_values)
411 , A_block_rowptr(block_graph.row_map)
412 , A_point_rowptr(point_graph.row_map)
413 , A_colind(block_graph.entries)
414 , blocksize_requested(blocksize_requested_)
415 , part2packrowidx0(interf.part2packrowidx0)
416 , part2rowidx0(interf.part2rowidx0)
417 , rowidx2part(interf.rowidx2part)
418 , partptr(interf.partptr)
419 , lclrow(interf.lclrow)
421 , A_x_offsets(amd.A_x_offsets)
422 , A_x_offsets_remote(amd.A_x_offsets_remote)
423 , is_dm2cm_active(dm2cm_.span() > 0)
424 , hasBlockCrsMatrix(hasBlockCrsMatrix_) {}
427 SerialDot(
const local_ordinal_type &blocksize,
428 const local_ordinal_type &lclRowID,
429 const local_ordinal_type &lclColID,
430 const local_ordinal_type &ii,
431 const ConstUnmanaged<local_ordinal_type_1d_view> colindsub_,
432 const impl_scalar_type *
const KOKKOS_RESTRICT xx,
433 impl_scalar_type *KOKKOS_RESTRICT yy)
const {
434 const size_type Aj_c = colindsub_(lclColID);
435 auto point_row_offset = A_point_rowptr(lclRowID * blocksize + ii) + Aj_c * blocksize;
436 impl_scalar_type val = 0;
437#if defined(KOKKOS_ENABLE_PRAGMA_IVDEP)
440#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
443 for (local_ordinal_type k1 = 0; k1 < blocksize; ++k1)
444 val += tpetra_values(point_row_offset + k1) * xx[k1];
449 SerialGemv(
const local_ordinal_type &blocksize,
450 const impl_scalar_type *
const KOKKOS_RESTRICT AA,
451 const impl_scalar_type *
const KOKKOS_RESTRICT xx,
452 impl_scalar_type *KOKKOS_RESTRICT yy)
const {
453 using tlb = BlockHelperDetails::TpetraLittleBlock<Tpetra::Impl::BlockCrsMatrixLittleBlockArrayLayout>;
454 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0) {
455 impl_scalar_type val = 0;
456#if defined(KOKKOS_ENABLE_PRAGMA_IVDEP)
459#if defined(KOKKOS_ENABLE_PRAGMA_UNROLL)
462 for (local_ordinal_type k1 = 0; k1 < blocksize; ++k1)
463 val += AA[tlb::getFlatIndex(k0, k1, blocksize)] * xx[k1];
468 template <
typename bbViewType,
typename yyViewType>
469 KOKKOS_INLINE_FUNCTION
void
470 VectorCopy(
const member_type &member,
471 const local_ordinal_type &blocksize,
472 const bbViewType &bb,
473 const yyViewType &yy)
const {
474 Kokkos::parallel_for(Kokkos::ThreadVectorRange(member, blocksize), [&](
const local_ordinal_type &k0) {
475 yy(k0) =
static_cast<typename yyViewType::const_value_type
>(bb(k0));
479 template <
typename xxViewType,
typename yyViewType>
480 KOKKOS_INLINE_FUNCTION
void
481 VectorDot(
const member_type &member,
482 const local_ordinal_type &blocksize,
483 const local_ordinal_type &lclRowID,
484 const local_ordinal_type &lclColID,
485 const local_ordinal_type &ii,
486 const ConstUnmanaged<local_ordinal_type_1d_view> colindsub_,
487 const xxViewType &xx,
488 const yyViewType &yy)
const {
489 const size_type Aj_c = colindsub_(lclColID);
490 auto point_row_offset = A_point_rowptr(lclRowID * blocksize + ii) + Aj_c * blocksize;
491 impl_scalar_type val = 0;
492 Kokkos::parallel_reduce(
493 Kokkos::ThreadVectorRange(member, blocksize),
494 [&](
const local_ordinal_type &k1, impl_scalar_type &update) {
495 update += tpetra_values(point_row_offset + k1) * xx(k1);
498 Kokkos::single(Kokkos::PerThread(member),
500 Kokkos::atomic_add(&yy(ii),
typename yyViewType::const_value_type(-val));
505 template <
typename AAViewType,
typename xxViewType,
typename yyViewType>
506 KOKKOS_INLINE_FUNCTION
void
507 VectorGemv(
const member_type &member,
508 const local_ordinal_type &blocksize,
509 const AAViewType &AA,
510 const xxViewType &xx,
511 const yyViewType &yy)
const {
512 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0) {
513 impl_scalar_type val = 0;
514 Kokkos::parallel_reduce(
515 Kokkos::ThreadVectorRange(member, blocksize),
516 [&](
const local_ordinal_type &k1, impl_scalar_type &update) {
517 update += AA(k0, k1) * xx(k1);
520 Kokkos::single(Kokkos::PerThread(member),
522 Kokkos::atomic_add(&yy(k0), -val);
529 KOKKOS_INLINE_FUNCTION
531 operator()(
const SeqTag &,
const local_ordinal_type &i)
const {
532 const local_ordinal_type blocksize = blocksize_requested;
533 const local_ordinal_type blocksize_square = blocksize * blocksize;
536 const Kokkos::pair<local_ordinal_type, local_ordinal_type> block_range(0, blocksize);
537 const local_ordinal_type num_vectors = y.extent(1);
538 const local_ordinal_type row = i * blocksize;
539 for (local_ordinal_type col = 0; col < num_vectors; ++col) {
541 impl_scalar_type *yy = &y(row, col);
542 const impl_scalar_type *
const bb = &b(row, col);
543 memcpy(yy, bb,
sizeof(impl_scalar_type) * blocksize);
546 const size_type A_k0 = A_block_rowptr[i];
547 for (size_type k = rowptr[i]; k < rowptr[i + 1]; ++k) {
548 const size_type j = A_k0 + colindsub[k];
549 const impl_scalar_type *
const xx = &x(A_colind[j] * blocksize, col);
550 if (hasBlockCrsMatrix) {
551 const impl_scalar_type *
const AA = &tpetra_values(j * blocksize_square);
552 SerialGemv(blocksize, AA, xx, yy);
554 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
555 SerialDot(blocksize, i, k, k0, colindsub, xx, yy);
561 KOKKOS_INLINE_FUNCTION
563 operator()(
const SeqTag &,
const member_type &member)
const {
565 const local_ordinal_type blocksize = blocksize_requested;
566 const local_ordinal_type blocksize_square = blocksize * blocksize;
568 const local_ordinal_type lr = member.league_rank();
569 const Kokkos::pair<local_ordinal_type, local_ordinal_type> block_range(0, blocksize);
570 const local_ordinal_type num_vectors = y.extent(1);
573 auto bb = Kokkos::subview(b, block_range, 0);
575 auto A_block_cst = ConstUnmanaged<tpetra_block_access_view_type>(tpetra_values.data(), blocksize, blocksize);
577 const local_ordinal_type row = lr * blocksize;
578 for (local_ordinal_type col = 0; col < num_vectors; ++col) {
580 auto yy = Kokkos::subview(y, Kokkos::make_pair(row, row + blocksize), col);
581 bb.assign_data(&b(row, col));
582 if (member.team_rank() == 0)
583 VectorCopy(member, blocksize, bb, yy);
584 member.team_barrier();
587 const size_type A_k0 = A_block_rowptr[lr];
589 if (hasBlockCrsMatrix) {
590 Kokkos::parallel_for(Kokkos::TeamThreadRange(member, rowptr[lr], rowptr[lr + 1]),
591 [&](
const local_ordinal_type &k) {
592 const size_type j = A_k0 + colindsub[k];
593 xx.assign_data(&x(A_colind[j] * blocksize, col));
594 A_block_cst.assign_data(&tpetra_values(j * blocksize_square));
595 VectorGemv(member, blocksize, A_block_cst, xx, yy);
598 Kokkos::parallel_for(Kokkos::TeamThreadRange(member, rowptr[lr], rowptr[lr + 1]),
599 [&](
const local_ordinal_type &k) {
600 const size_type j = A_k0 + colindsub[k];
601 xx.assign_data(&x(A_colind[j] * blocksize, col));
603 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
604 VectorDot(member, blocksize, lr, k, k0, colindsub, xx, yy);
616 template <
int B,
bool async,
bool overlap,
bool haveBlockMatrix>
618 static_assert(!(async && overlap),
619 "ComputeResidualVector: async && overlap is not a valid configuration for GeneralTag");
625 template <
int P,
int B,
bool haveBlockMatrix>
626 using OverlapTag = GeneralTag<B, false, P != 0, haveBlockMatrix>;
628 template <
int B,
bool haveBlockMatrix>
629 using AsyncTag = GeneralTag<B, true, false, haveBlockMatrix>;
632 template <
int B,
bool async,
bool overlap,
bool haveBlockMatrix>
633 KOKKOS_INLINE_FUNCTION
void
634 operator()(
const GeneralTag<B, async, overlap, haveBlockMatrix> &,
const local_ordinal_type &rowidx)
const {
635 const local_ordinal_type blocksize = (B == 0 ? blocksize_requested : B);
638 const local_ordinal_type partidx = rowidx2part(rowidx);
639 const local_ordinal_type pri = part2packrowidx0(partidx) + (rowidx - partptr(partidx));
640 const local_ordinal_type v = partidx % vector_length;
642 const local_ordinal_type num_vectors = y_packed.extent(2);
643 const local_ordinal_type num_local_rows = lclrow.extent(0);
646 impl_scalar_type yy[B == 0 ? max_blocksize : B] = {};
648 const local_ordinal_type lr = lclrow(rowidx);
650 auto colindsub_used = overlap ? colindsub_remote : colindsub;
651 auto rowptr_used = overlap ? rowptr_remote : rowptr;
653 for (local_ordinal_type col = 0; col < num_vectors; ++col) {
654 if constexpr (overlap) {
656 memset((
void *)yy, 0,
sizeof(impl_scalar_type) * blocksize);
659 const local_ordinal_type row = lr * blocksize;
660 memcpy(yy, &b(row, col),
sizeof(impl_scalar_type) * blocksize);
664 const size_type A_k0 = A_block_rowptr[lr];
665 for (size_type k = rowptr_used[lr]; k < rowptr_used[lr + 1]; ++k) {
666 const size_type j = A_k0 + colindsub_used[k];
667 const local_ordinal_type A_colind_at_j = A_colind[j];
668 if constexpr (haveBlockMatrix) {
669 const local_ordinal_type blocksize_square = blocksize * blocksize;
670 const impl_scalar_type *
const AA = &tpetra_values(j * blocksize_square);
671 if ((!async && !overlap) || (async && A_colind_at_j < num_local_rows)) {
672 const auto loc = is_dm2cm_active ? dm2cm[A_colind_at_j] : A_colind_at_j;
673 const impl_scalar_type *
const xx = &x(loc * blocksize, col);
674 SerialGemv(blocksize, AA, xx, yy);
676 const auto loc = A_colind_at_j - num_local_rows;
677 const impl_scalar_type *
const xx_remote = &x_remote(loc * blocksize, col);
678 SerialGemv(blocksize, AA, xx_remote, yy);
681 if ((!async && !overlap) || (async && A_colind_at_j < num_local_rows)) {
682 const auto loc = is_dm2cm_active ? dm2cm[A_colind_at_j] : A_colind_at_j;
683 const impl_scalar_type *
const xx = &x(loc * blocksize, col);
684 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
685 SerialDot(blocksize, lr, k, k0, colindsub_used, xx, yy);
687 const auto loc = A_colind_at_j - num_local_rows;
688 const impl_scalar_type *
const xx_remote = &x_remote(loc * blocksize, col);
689 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
690 SerialDot(blocksize, lr, k, k0, colindsub_used, xx_remote, yy);
695 if constexpr (overlap) {
696 for (local_ordinal_type k = 0; k < blocksize; ++k)
697 y_packed(pri, k, col)[v] += yy[k];
699 for (local_ordinal_type k = 0; k < blocksize; ++k)
700 y_packed(pri, k, col)[v] = yy[k];
706 template <
int B,
bool async,
bool overlap>
707 KOKKOS_INLINE_FUNCTION
void
708 operator()(
const GeneralTag<B, async, overlap, true> &,
const member_type &member)
const {
709 const local_ordinal_type blocksize = (B == 0 ? blocksize_requested : B);
712 const local_ordinal_type rowidx = member.league_rank();
713 const local_ordinal_type partidx = rowidx2part(rowidx);
714 const local_ordinal_type pri = part2packrowidx0(partidx) + (rowidx - partptr(partidx));
715 const local_ordinal_type v = partidx % vector_length;
717 const Kokkos::pair<local_ordinal_type, local_ordinal_type> block_range(0, blocksize);
718 const local_ordinal_type num_vectors = y_packed_scalar.extent(2);
719 const local_ordinal_type num_local_rows = lclrow.extent(0);
722 auto bb = Kokkos::subview(b, block_range, 0);
724 auto yy = Kokkos::subview(y_packed_scalar, 0, block_range, 0, 0);
725 auto A_block_cst = ConstUnmanaged<tpetra_block_access_view_type>(tpetra_values.data(), blocksize, blocksize);
728 impl_scalar_type *local_Ax =
reinterpret_cast<impl_scalar_type *
>(member.team_scratch(0).get_shmem(blocksize *
sizeof(impl_scalar_type)));
729 impl_scalar_type *local_x =
reinterpret_cast<impl_scalar_type *
>(member.thread_scratch(0).get_shmem(blocksize *
sizeof(impl_scalar_type)));
731 const local_ordinal_type lr = lclrow(rowidx);
732 const local_ordinal_type row = lr * blocksize;
733 for (local_ordinal_type col = 0; col < num_vectors; ++col) {
735 member.team_barrier();
738 Kokkos::parallel_for(Kokkos::TeamVectorRange(member, blocksize), [&](
const local_ordinal_type &i) {
741 member.team_barrier();
744 if constexpr (!overlap) {
745 numEntries = A_x_offsets.extent(2);
747 numEntries = A_x_offsets_remote.extent(2);
750 Kokkos::parallel_for(Kokkos::TeamThreadRange(member, 0, numEntries),
752 int64_t A_offset = overlap ? A_x_offsets_remote(rowidx, 0, k) : A_x_offsets(rowidx, 0, k);
753 int64_t x_offset = overlap ? A_x_offsets_remote(rowidx, 1, k) : A_x_offsets(rowidx, 1, k);
754 if (A_offset != KokkosKernels::ArithTraits<int64_t>::min()) {
755 A_block_cst.assign_data(tpetra_values.data() + A_offset);
757 if constexpr (async) {
758 size_type remote_cutoff = blocksize * num_local_rows;
759 if (x_offset >= remote_cutoff)
760 xx.assign_data(&x_remote(x_offset - remote_cutoff, col));
762 xx.assign_data(&x(x_offset, col));
764 if constexpr (!overlap) {
765 xx.assign_data(&x(x_offset, col));
767 xx.assign_data(&x_remote(x_offset, col));
771 Kokkos::parallel_for(Kokkos::ThreadVectorRange(member, blocksize), [&](
const local_ordinal_type &i) {
776 Kokkos::parallel_for(
777 Kokkos::ThreadVectorRange(member, blocksize),
778 [&](
const local_ordinal_type &k0) {
779 impl_scalar_type val = 0;
780 for (
int k1 = 0; k1 < blocksize; k1++)
781 val += A_block_cst(k0, k1) * local_x[k1];
782 Kokkos::atomic_add(local_Ax + k0, val);
786 member.team_barrier();
788 yy.assign_data(&y_packed_scalar(pri, 0, col, v));
789 bb.assign_data(&b(row, col));
790 Kokkos::parallel_for(Kokkos::TeamVectorRange(member, blocksize), [&](
const local_ordinal_type &i) {
792 yy(i) = bb(i) - local_Ax[i];
794 yy(i) -= local_Ax[i];
800 template <
int B,
bool async,
bool overlap>
801 KOKKOS_INLINE_FUNCTION
void
802 operator()(
const GeneralTag<B, async, overlap, false> &,
const member_type &member)
const {
803 const local_ordinal_type blocksize = (B == 0 ? blocksize_requested : B);
806 const local_ordinal_type rowidx = member.league_rank();
807 const local_ordinal_type partidx = rowidx2part(rowidx);
808 const local_ordinal_type pri = part2packrowidx0(partidx) + (rowidx - partptr(partidx));
809 const local_ordinal_type v = partidx % vector_length;
811 const Kokkos::pair<local_ordinal_type, local_ordinal_type> block_range(0, blocksize);
812 const local_ordinal_type num_vectors = y_packed_scalar.extent(2);
813 const local_ordinal_type num_local_rows = lclrow.extent(0);
816 auto bb = Kokkos::subview(b, block_range, 0);
819 auto yy = Kokkos::subview(y_packed_scalar, 0, block_range, 0, 0);
820 auto A_block_cst = ConstUnmanaged<tpetra_block_access_view_type>(tpetra_values.data(), blocksize, blocksize);
821 auto colindsub_used = overlap ? colindsub_remote : colindsub;
822 auto rowptr_used = overlap ? rowptr_remote : rowptr;
824 const local_ordinal_type lr = lclrow(rowidx);
825 const local_ordinal_type row = lr * blocksize;
826 for (local_ordinal_type col = 0; col < num_vectors; ++col) {
827 yy.assign_data(&y_packed_scalar(pri, 0, col, v));
830 bb.assign_data(&b(row, col));
831 if (member.team_rank() == 0)
832 VectorCopy(member, blocksize, bb, yy);
833 member.team_barrier();
837 const size_type A_k0 = A_block_rowptr[lr];
838 Kokkos::parallel_for(Kokkos::TeamThreadRange(member, rowptr_used[lr], rowptr_used[lr + 1]),
839 [&](
const local_ordinal_type &k) {
840 const size_type j = A_k0 + colindsub_used[k];
841 const local_ordinal_type A_colind_at_j = A_colind[j];
842 if ((async && A_colind_at_j < num_local_rows) || (!async && !overlap)) {
843 const auto loc = is_dm2cm_active ? dm2cm[A_colind_at_j] : A_colind_at_j;
844 xx.assign_data(&x(loc * blocksize, col));
845 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
846 VectorDot(member, blocksize, lr, k, k0, colindsub_used, xx, yy);
848 const auto loc = A_colind_at_j - num_local_rows;
849 xx_remote.assign_data(&x_remote(loc * blocksize, col));
850 for (local_ordinal_type k0 = 0; k0 < blocksize; ++k0)
851 VectorDot(member, blocksize, lr, k, k0, colindsub_used, xx_remote, yy);
858 template <
typename MultiVectorLocalViewTypeY,
859 typename MultiVectorLocalViewTypeB,
860 typename MultiVectorLocalViewTypeX>
861 void run(
const MultiVectorLocalViewTypeY &y_,
862 const MultiVectorLocalViewTypeB &b_,
863 const MultiVectorLocalViewTypeX &x_) {
864 IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN;
865 IFPACK2_BLOCKHELPER_TIMER_WITH_FENCE(
"BlockTriDi::ComputeResidual::<SeqTag>", ComputeResidual0, execution_space);
870 if constexpr (is_device<execution_space>::value) {
871 const local_ordinal_type blocksize = blocksize_requested;
872 const local_ordinal_type team_size = 8;
873 const local_ordinal_type vector_size = ComputeResidualVectorRecommendedVectorSize<execution_space>(blocksize, team_size);
874 const Kokkos::TeamPolicy<execution_space, SeqTag> policy(rowptr.extent(0) - 1, team_size, vector_size);
875 Kokkos::parallel_for(
"ComputeResidual::TeamPolicy::run<SeqTag>", policy, *
this);
877 const Kokkos::RangePolicy<execution_space, SeqTag> policy(0, rowptr.extent(0) - 1);
878 Kokkos::parallel_for(
"ComputeResidual::RangePolicy::run<SeqTag>", policy, *
this);
880 IFPACK2_BLOCKHELPER_PROFILER_REGION_END;
881 IFPACK2_BLOCKHELPER_TIMER_FENCE(execution_space)
885 template <
typename MultiVectorLocalViewTypeB,
886 typename MultiVectorLocalViewTypeX,
887 typename MultiVectorLocalViewTypeX_Remote>
888 void run(
const vector_type_3d_view &y_packed_,
889 const MultiVectorLocalViewTypeB &b_,
890 const MultiVectorLocalViewTypeX &x_,
891 const MultiVectorLocalViewTypeX_Remote &x_remote_) {
892 IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN;
893 IFPACK2_BLOCKHELPER_TIMER_WITH_FENCE(
"BlockTriDi::ComputeResidual::<AsyncTag>", ComputeResidual0, execution_space);
897 x_remote = x_remote_;
898 if constexpr (is_device<execution_space>::value) {
899 y_packed_scalar = btdm_scalar_type_4d_view((btdm_scalar_type *)y_packed_.data(),
905 y_packed = y_packed_;
908 if constexpr (is_device<execution_space>::value) {
909 const local_ordinal_type blocksize = blocksize_requested;
914#define BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(B) \
916 if (this->hasBlockCrsMatrix) { \
917 const local_ordinal_type team_size = 8; \
918 const local_ordinal_type vector_size = 8; \
919 const size_t shmem_team_size = blocksize * sizeof(btdm_scalar_type); \
920 const size_t shmem_thread_size = blocksize * sizeof(btdm_scalar_type); \
921 Kokkos::TeamPolicy<execution_space, AsyncTag<B, true>> \
922 policy(rowidx2part.extent(0), team_size, vector_size); \
923 policy.set_scratch_size(0, Kokkos::PerTeam(shmem_team_size), Kokkos::PerThread(shmem_thread_size)); \
924 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<AsyncTag>", \
927 const local_ordinal_type team_size = 8; \
928 const local_ordinal_type vector_size = ComputeResidualVectorRecommendedVectorSize<execution_space>(blocksize, team_size); \
929 const Kokkos::TeamPolicy<execution_space, AsyncTag<B, false>> \
930 policy(rowidx2part.extent(0), team_size, vector_size); \
931 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<AsyncTag>", \
936 switch (blocksize_requested) {
937 case 3: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(3);
938 case 5: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(5);
939 case 7: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(7);
940 case 9: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(9);
941 case 10: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(10);
942 case 11: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(11);
943 case 16: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(16);
944 case 17: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(17);
945 case 18: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(18);
946 default: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(0);
948#undef BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL
950#define BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(B) \
952 if (this->hasBlockCrsMatrix) { \
953 const Kokkos::RangePolicy<execution_space, AsyncTag<B, true>> policy(0, rowidx2part.extent(0)); \
954 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<AsyncTag>", \
957 const Kokkos::RangePolicy<execution_space, AsyncTag<B, false>> policy(0, rowidx2part.extent(0)); \
958 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<AsyncTag>", \
964 switch (blocksize_requested) {
965 case 3: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(3);
966 case 5: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(5);
967 case 7: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(7);
968 case 9: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(9);
969 case 10: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(10);
970 case 11: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(11);
971 case 16: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(16);
972 case 17: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(17);
973 case 18: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(18);
974 default: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(0);
976#undef BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL
978 IFPACK2_BLOCKHELPER_PROFILER_REGION_END;
979 IFPACK2_BLOCKHELPER_TIMER_FENCE(execution_space)
983 template <
typename MultiVectorLocalViewTypeB,
984 typename MultiVectorLocalViewTypeX,
985 typename MultiVectorLocalViewTypeX_Remote>
986 void run(
const vector_type_3d_view &y_packed_,
987 const MultiVectorLocalViewTypeB &b_,
988 const MultiVectorLocalViewTypeX &x_,
989 const MultiVectorLocalViewTypeX_Remote &x_remote_,
990 const bool compute_owned) {
991 IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN;
992 IFPACK2_BLOCKHELPER_TIMER_WITH_FENCE(
"BlockTriDi::ComputeResidual::<OverlapTag>", ComputeResidual0, execution_space);
996 x_remote = x_remote_;
997 if constexpr (is_device<execution_space>::value) {
998 y_packed_scalar = btdm_scalar_type_4d_view((btdm_scalar_type *)y_packed_.data(),
1000 y_packed_.extent(1),
1001 y_packed_.extent(2),
1004 y_packed = y_packed_;
1007 if constexpr (is_device<execution_space>::value) {
1008 const local_ordinal_type blocksize = blocksize_requested;
1013#define BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(B) \
1014 if (this->hasBlockCrsMatrix) { \
1015 const local_ordinal_type team_size = 8; \
1016 const local_ordinal_type vector_size = 8; \
1017 const size_t shmem_team_size = blocksize * sizeof(btdm_scalar_type); \
1018 const size_t shmem_thread_size = blocksize * sizeof(btdm_scalar_type); \
1019 if (compute_owned) { \
1020 Kokkos::TeamPolicy<execution_space, OverlapTag<0, B, true>> \
1021 policy(rowidx2part.extent(0), team_size, vector_size); \
1022 policy.set_scratch_size(0, Kokkos::PerTeam(shmem_team_size), Kokkos::PerThread(shmem_thread_size)); \
1023 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<OverlapTag<0> >", policy, *this); \
1025 Kokkos::TeamPolicy<execution_space, OverlapTag<1, B, true>> \
1026 policy(rowidx2part.extent(0), team_size, vector_size); \
1027 policy.set_scratch_size(0, Kokkos::PerTeam(shmem_team_size), Kokkos::PerThread(shmem_thread_size)); \
1028 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<OverlapTag<1> >", policy, *this); \
1031 const local_ordinal_type team_size = 8; \
1032 const local_ordinal_type vector_size = ComputeResidualVectorRecommendedVectorSize<execution_space>(blocksize, team_size); \
1033 if (compute_owned) { \
1034 const Kokkos::TeamPolicy<execution_space, OverlapTag<0, B, false>> \
1035 policy(rowidx2part.extent(0), team_size, vector_size); \
1036 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<OverlapTag<0> >", policy, *this); \
1038 const Kokkos::TeamPolicy<execution_space, OverlapTag<1, B, false>> \
1039 policy(rowidx2part.extent(0), team_size, vector_size); \
1040 Kokkos::parallel_for("ComputeResidual::TeamPolicy::run<OverlapTag<1> >", policy, *this); \
1044 switch (blocksize_requested) {
1045 case 3: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(3);
1046 case 5: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(5);
1047 case 7: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(7);
1048 case 9: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(9);
1049 case 10: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(10);
1050 case 11: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(11);
1051 case 16: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(16);
1052 case 17: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(17);
1053 case 18: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(18);
1054 default: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(0);
1056#undef BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL
1058#define BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(B) \
1059 if (this->hasBlockCrsMatrix) { \
1060 if (compute_owned) { \
1061 const Kokkos::RangePolicy<execution_space, OverlapTag<0, B, true>> \
1062 policy(0, rowidx2part.extent(0)); \
1063 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<OverlapTag<0> >", policy, *this); \
1065 const Kokkos::RangePolicy<execution_space, OverlapTag<1, B, true>> \
1066 policy(0, rowidx2part.extent(0)); \
1067 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<OverlapTag<1> >", policy, *this); \
1070 if (compute_owned) { \
1071 const Kokkos::RangePolicy<execution_space, OverlapTag<0, B, false>> \
1072 policy(0, rowidx2part.extent(0)); \
1073 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<OverlapTag<0> >", policy, *this); \
1075 const Kokkos::RangePolicy<execution_space, OverlapTag<1, B, false>> \
1076 policy(0, rowidx2part.extent(0)); \
1077 Kokkos::parallel_for("ComputeResidual::RangePolicy::run<OverlapTag<1> >", policy, *this); \
1082 switch (blocksize_requested) {
1083 case 3: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(3);
1084 case 5: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(5);
1085 case 7: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(7);
1086 case 9: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(9);
1087 case 10: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(10);
1088 case 11: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(11);
1089 case 16: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(16);
1090 case 17: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(17);
1091 case 18: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(18);
1092 default: BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL(0);
1094#undef BLOCKTRIDICONTAINER_DETAILS_COMPUTERESIDUAL
1096 IFPACK2_BLOCKHELPER_PROFILER_REGION_END;
1097 IFPACK2_BLOCKHELPER_TIMER_FENCE(execution_space)
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Definition Ifpack2_BlockComputeResidualVector.hpp:23
Definition Ifpack2_BlockHelper.hpp:270
size_t size_type
Definition Ifpack2_BlockHelper.hpp:274
node_type::device_type node_device_type
Definition Ifpack2_BlockHelper.hpp:297
KokkosKernels::ArithTraits< scalar_type >::val_type impl_scalar_type
Definition Ifpack2_BlockHelper.hpp:283
Kokkos::View< size_type *, device_type > size_type_1d_view
Definition Ifpack2_BlockHelper.hpp:346