Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_BlockHelper.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_BLOCKHELPER_IMPL_HPP
11#define IFPACK2_BLOCKHELPER_IMPL_HPP
12
13#include "Ifpack2_BlockHelper_Timers.hpp"
14
15namespace Ifpack2 {
16
17namespace BlockHelperDetails {
18
19namespace KB = KokkosBatched;
20
24using do_not_initialize_tag = Kokkos::ViewAllocateWithoutInitializing;
25
26template <typename MemoryTraitsType, Kokkos::MemoryTraitsFlags flag>
27using MemoryTraits = Kokkos::MemoryTraits<MemoryTraitsType::is_unmanaged |
28 MemoryTraitsType::is_random_access |
29 flag>;
30
31template <typename ViewType>
32using Unmanaged = Kokkos::View<typename ViewType::data_type,
33 typename ViewType::array_layout,
34 typename ViewType::device_type,
35 MemoryTraits<typename ViewType::memory_traits, Kokkos::Unmanaged> >;
36template <typename ViewType>
37using Atomic = Kokkos::View<typename ViewType::data_type,
38 typename ViewType::array_layout,
39 typename ViewType::device_type,
40 MemoryTraits<typename ViewType::memory_traits, Kokkos::Atomic> >;
41template <typename ViewType>
42using Const = Kokkos::View<typename ViewType::const_data_type,
43 typename ViewType::array_layout,
44 typename ViewType::device_type,
45 typename ViewType::memory_traits>;
46template <typename ViewType>
47using ConstUnmanaged = Const<Unmanaged<ViewType> >;
48
49template <typename ViewType>
50using AtomicUnmanaged = Atomic<Unmanaged<ViewType> >;
51
52template <typename ViewType>
53using Unmanaged = Kokkos::View<typename ViewType::data_type,
54 typename ViewType::array_layout,
55 typename ViewType::device_type,
56 MemoryTraits<typename ViewType::memory_traits, Kokkos::Unmanaged> >;
57
58template <typename ViewType>
59using Scratch = Kokkos::View<typename ViewType::data_type,
60 typename ViewType::array_layout,
61 typename ViewType::execution_space::scratch_memory_space,
62 MemoryTraits<typename ViewType::memory_traits, Kokkos::Unmanaged> >;
63
67template <typename LayoutType>
69template <>
70struct TpetraLittleBlock<Kokkos::LayoutLeft> {
71 template <typename T>
72 KOKKOS_INLINE_FUNCTION static T getFlatIndex(const T i, const T j, const T blksize) { return i + j * blksize; }
73};
74template <>
75struct TpetraLittleBlock<Kokkos::LayoutRight> {
76 template <typename T>
77 KOKKOS_INLINE_FUNCTION static T getFlatIndex(const T i, const T j, const T blksize) { return i * blksize + j; }
78};
79
83template <typename T>
84struct BlockTridiagScalarType { typedef T type; };
85#if defined(IFPACK2_BLOCKHELPER_USE_SMALL_SCALAR_FOR_BLOCKTRIDIAG)
86template <>
87struct BlockTridiagScalarType<double> { typedef float type; };
88// template<> struct SmallScalarType<Kokkos::complex<double> > { typedef Kokkos::complex<float> type; };
89#endif
90
94template <typename T>
95struct is_cuda {
96 enum : bool { value = false };
97};
98#if defined(KOKKOS_ENABLE_CUDA)
99template <>
100struct is_cuda<Kokkos::Cuda> {
101 enum : bool { value = true };
102};
103#endif
104
108template <typename T>
109struct is_hip {
110 enum : bool { value = false };
111};
112#if defined(KOKKOS_ENABLE_HIP)
113template <>
114struct is_hip<Kokkos::HIP> {
115 enum : bool { value = true };
116};
117#endif
118
122template <typename T>
123struct is_sycl {
124 enum : bool { value = false };
125};
126#if defined(KOKKOS_ENABLE_SYCL)
127template <>
128struct is_sycl<Kokkos::Experimental::SYCL> {
129 enum : bool { value = true };
130};
131#endif
132
133template <typename T>
134struct is_device {
135 enum : bool { value = is_cuda<T>::value || is_hip<T>::value || is_sycl<T>::value };
136};
137
141template <typename T>
143 static void createInstance(T &exec_instance) {
144 exec_instance = T();
145 }
146#if defined(KOKKOS_ENABLE_CUDA)
147 static void createInstance(const cudaStream_t &s, T &exec_instance) {
148 exec_instance = T();
149 }
150#endif
151};
152
153#if defined(KOKKOS_ENABLE_CUDA)
154template <>
155struct ExecutionSpaceFactory<Kokkos::Cuda> {
156 static void createInstance(Kokkos::Cuda &exec_instance) {
157 exec_instance = Kokkos::Cuda();
158 }
159 static void createInstance(const cudaStream_t &s, Kokkos::Cuda &exec_instance) {
160 exec_instance = Kokkos::Cuda(s);
161 }
162};
163#endif
164
165#if defined(KOKKOS_ENABLE_HIP)
166template <>
167struct ExecutionSpaceFactory<Kokkos::HIP> {
168 static void createInstance(Kokkos::HIP &exec_instance) {
169 exec_instance = Kokkos::HIP();
170 }
171};
172#endif
173
174#if defined(KOKKOS_ENABLE_SYCL)
175template <>
176struct ExecutionSpaceFactory<Kokkos::Experimental::SYCL> {
177 static void createInstance(Kokkos::Experimental::SYCL &exec_instance) {
178 exec_instance = Kokkos::Experimental::SYCL();
179 }
180};
181#endif
182
183#if defined(KOKKOS_ENABLE_CUDA) && defined(IFPACK2_BLOCKHELPER_ENABLE_PROFILE)
184#define IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN \
185 KOKKOS_IMPL_CUDA_SAFE_CALL(cudaProfilerStart());
186
187#define IFPACK2_BLOCKHELPER_PROFILER_REGION_END \
188 { KOKKOS_IMPL_CUDA_SAFE_CALL(cudaProfilerStop()); }
189#else
191#define IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN
192#define IFPACK2_BLOCKHELPER_PROFILER_REGION_END
193#endif
194
198template <typename CommPtrType>
199std::string get_msg_prefix(const CommPtrType &comm) {
200 const auto rank = comm->getRank();
201 const auto nranks = comm->getSize();
202 std::stringstream ss;
203 ss << "Rank " << rank << " of " << nranks << ": ";
204 return ss.str();
205}
206
210template <typename T, int N>
212 T v[N];
213 KOKKOS_INLINE_FUNCTION
215 for (int i = 0; i < N; ++i)
216 this->v[i] = 0;
217 }
218 KOKKOS_INLINE_FUNCTION
220 for (int i = 0; i < N; ++i)
221 this->v[i] = b.v[i];
222 }
223};
224template <typename T, int N>
225static KOKKOS_INLINE_FUNCTION void
226operator+=(ArrayValueType<T, N> &a,
227 const ArrayValueType<T, N> &b) {
228 for (int i = 0; i < N; ++i)
229 a.v[i] += b.v[i];
230}
231
235template <typename T, int N, typename ExecSpace>
237 typedef SumReducer reducer;
239 typedef Kokkos::View<value_type, ExecSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > result_view_type;
240 value_type *value;
241
242 KOKKOS_INLINE_FUNCTION
244 : value(&val) {}
245
246 KOKKOS_INLINE_FUNCTION
247 void join(value_type &dst, value_type const &src) const {
248 for (int i = 0; i < N; ++i)
249 dst.v[i] += src.v[i];
250 }
251 KOKKOS_INLINE_FUNCTION
252 void init(value_type &val) const {
253 for (int i = 0; i < N; ++i)
254 val.v[i] = 0;
255 }
256 KOKKOS_INLINE_FUNCTION
257 value_type &reference() {
258 return *value;
259 }
260 KOKKOS_INLINE_FUNCTION
261 result_view_type view() const {
262 return result_view_type(value);
263 }
264};
265
269template <typename MatrixType>
270struct ImplType {
274 typedef size_t size_type;
275 typedef typename MatrixType::scalar_type scalar_type;
276 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
277 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
278 typedef typename MatrixType::node_type node_type;
279
283 typedef typename KokkosKernels::ArithTraits<scalar_type>::val_type impl_scalar_type;
284 typedef typename KokkosKernels::ArithTraits<impl_scalar_type>::mag_type magnitude_type;
285
286 typedef typename BlockTridiagScalarType<impl_scalar_type>::type btdm_scalar_type;
287 typedef typename KokkosKernels::ArithTraits<btdm_scalar_type>::mag_type btdm_magnitude_type;
288
292 typedef Kokkos::DefaultHostExecutionSpace host_execution_space;
293
297 typedef typename node_type::device_type node_device_type;
298 typedef typename node_device_type::execution_space node_execution_space;
299 typedef typename node_device_type::memory_space node_memory_space;
300
301#if defined(KOKKOS_ENABLE_CUDA) && defined(IFPACK2_BLOCKHELPER_USE_CUDA_SPACE)
303 typedef node_execution_space execution_space;
304 typedef typename std::conditional<std::is_same<node_memory_space, Kokkos::CudaUVMSpace>::value,
305 Kokkos::CudaSpace,
306 node_memory_space>::type memory_space;
307 typedef Kokkos::Device<execution_space, memory_space> device_type;
308#else
309 typedef node_device_type device_type;
310 typedef node_execution_space execution_space;
311 typedef node_memory_space memory_space;
312#endif
313
314 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> tpetra_multivector_type;
315 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> tpetra_map_type;
316 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> tpetra_import_type;
317 typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> tpetra_row_matrix_type;
318 typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> tpetra_crs_matrix_type;
319 typedef Tpetra::CrsGraph<local_ordinal_type, global_ordinal_type, node_type> tpetra_crs_graph_type;
320 typedef Tpetra::BlockCrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> tpetra_block_crs_matrix_type;
321 typedef typename tpetra_block_crs_matrix_type::little_block_type tpetra_block_access_view_type;
322 typedef Tpetra::BlockMultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> tpetra_block_multivector_type;
323 typedef typename tpetra_block_crs_matrix_type::crs_graph_type::local_graph_device_type local_crs_graph_type;
324
328 template <typename T, int l>
329 using Vector = KB::Vector<T, l>;
330 template <typename T>
331 using SIMD = KB::SIMD<T>;
332 template <typename T, typename M>
333 using DefaultVectorLength = KB::DefaultVectorLength<T, M>;
334 template <typename T, typename M>
335 using DefaultInternalVectorLength = KB::DefaultInternalVectorLength<T, M>;
336
337 static constexpr int vector_length = DefaultVectorLength<btdm_scalar_type, memory_space>::value;
338 static constexpr int internal_vector_length = DefaultInternalVectorLength<btdm_scalar_type, memory_space>::value;
339 static constexpr int half_vector_length = (vector_length > 1) ? (vector_length / 2) : 1;
340 typedef Vector<SIMD<btdm_scalar_type>, vector_length> vector_type;
341 typedef Vector<SIMD<btdm_scalar_type>, internal_vector_length> internal_vector_type;
342
346 typedef Kokkos::View<size_type *, device_type> size_type_1d_view;
347 typedef Kokkos::View<size_type **, device_type> size_type_2d_view;
348 typedef Kokkos::View<int64_t ***, Kokkos::LayoutRight, device_type> i64_3d_view;
349 typedef Kokkos::View<local_ordinal_type *, device_type> local_ordinal_type_1d_view;
350 typedef Kokkos::View<local_ordinal_type **, device_type> local_ordinal_type_2d_view;
351 // tpetra block crs values
352 typedef Kokkos::View<impl_scalar_type *, device_type> impl_scalar_type_1d_view;
353 typedef Kokkos::View<impl_scalar_type *, node_device_type> impl_scalar_type_1d_view_tpetra;
354
355 // tpetra multivector values (layout left): may need to change the typename more explicitly
356 typedef Kokkos::View<impl_scalar_type **, Kokkos::LayoutLeft, device_type> impl_scalar_type_2d_view;
357 typedef Kokkos::View<impl_scalar_type **, Kokkos::LayoutLeft, node_device_type> impl_scalar_type_2d_view_tpetra;
358 typedef Kokkos::View<const impl_scalar_type **, Kokkos::LayoutLeft, node_device_type> const_impl_scalar_type_2d_view_tpetra;
359
360 // packed data always use layout right
361 typedef Kokkos::View<vector_type *, device_type> vector_type_1d_view;
362 typedef Kokkos::View<vector_type ***, Kokkos::LayoutRight, device_type> vector_type_3d_view;
363 typedef Kokkos::View<vector_type ****, Kokkos::LayoutRight, device_type> vector_type_4d_view;
364 typedef Kokkos::View<internal_vector_type ***, Kokkos::LayoutRight, device_type> internal_vector_type_3d_view;
365 typedef Kokkos::View<internal_vector_type ****, Kokkos::LayoutRight, device_type> internal_vector_type_4d_view;
366 typedef Kokkos::View<internal_vector_type *****, Kokkos::LayoutRight, device_type> internal_vector_type_5d_view;
367 typedef Kokkos::View<btdm_scalar_type **, Kokkos::LayoutRight, device_type> btdm_scalar_type_2d_view;
368 typedef Kokkos::View<btdm_scalar_type ***, Kokkos::LayoutRight, device_type> btdm_scalar_type_3d_view;
369 typedef Kokkos::View<btdm_scalar_type ****, Kokkos::LayoutRight, device_type> btdm_scalar_type_4d_view;
370 typedef Kokkos::View<btdm_scalar_type *****, Kokkos::LayoutRight, device_type> btdm_scalar_type_5d_view;
371};
372
376template <typename MatrixType>
378 public:
380 using host_execution_space = typename impl_type::host_execution_space;
381 using magnitude_type = typename impl_type::magnitude_type;
382
383 private:
384 bool collective_;
385 int sweep_step_, sweep_step_upper_bound_;
386#ifdef HAVE_IFPACK2_MPI
387 MPI_Request mpi_request_;
388 MPI_Comm comm_;
389#endif
390 magnitude_type work_[3];
391
392 public:
393 NormManager() = default;
394 NormManager(const NormManager &b) = default;
395 NormManager(const Teuchos::RCP<const Teuchos::Comm<int> > &comm) {
396 sweep_step_ = 1;
397 sweep_step_upper_bound_ = 1;
398 collective_ = comm->getSize() > 1;
399 if (collective_) {
400#ifdef HAVE_IFPACK2_MPI
401 const auto mpi_comm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
402 TEUCHOS_ASSERT(!mpi_comm.is_null());
403 comm_ = *mpi_comm->getRawMpiComm();
404#endif
405 }
406 const magnitude_type zero(0), minus_one(-1);
407 work_[0] = zero;
408 work_[1] = zero;
409 work_[2] = minus_one;
410 }
411
412 // Check the norm every sweep_step sweeps.
413 void setCheckFrequency(const int sweep_step) {
414 TEUCHOS_TEST_FOR_EXCEPT_MSG(sweep_step < 1, "sweep step must be >= 1");
415 sweep_step_upper_bound_ = sweep_step;
416 sweep_step_ = 1;
417 }
418
419 // Get the buffer into which to store rank-local squared norms.
420 magnitude_type *getBuffer() { return &work_[0]; }
421
422 // Call MPI_Iallreduce to find the global squared norms.
423 void ireduce(const int sweep, const bool force = false) {
424 if (!force && sweep % sweep_step_) return;
425
426 IFPACK2_BLOCKHELPER_TIMER("BlockTriDi::NormManager::Ireduce", Ireduce);
427
428 work_[1] = work_[0];
429#ifdef HAVE_IFPACK2_MPI
430 auto send_data = &work_[1];
431 auto recv_data = &work_[0];
432 if (collective_) {
433#if defined(IFPACK2_BLOCKTRIDICONTAINER_USE_MPI_3)
434 MPI_Iallreduce(send_data, recv_data, 1,
435 Teuchos::Details::MpiTypeTraits<magnitude_type>::getType(),
436 MPI_SUM, comm_, &mpi_request_);
437#else
438 MPI_Allreduce(send_data, recv_data, 1,
439 Teuchos::Details::MpiTypeTraits<magnitude_type>::getType(),
440 MPI_SUM, comm_);
441#endif
442 }
443#endif
444 }
445
446 // Check if the norm-based termination criterion is met. tol2 is the
447 // tolerance squared. Sweep is the sweep index. If not every iteration is
448 // being checked, this function immediately returns false. If a check must
449 // be done at this iteration, it waits for the reduction triggered by
450 // ireduce to complete, then checks the global norm against the tolerance.
451 bool checkDone(const int sweep, const magnitude_type tol2, const bool force = false) {
452 // early return
453 if (sweep <= 0) return false;
454
455 IFPACK2_BLOCKHELPER_TIMER("BlockTriDi::NormManager::CheckDone", CheckDone);
456
457 TEUCHOS_ASSERT(sweep >= 1);
458 if (!force && (sweep - 1) % sweep_step_) return false;
459 if (collective_) {
460#ifdef HAVE_IFPACK2_MPI
461#if defined(IFPACK2_BLOCKTRIDICONTAINER_USE_MPI_3)
462 MPI_Wait(&mpi_request_, MPI_STATUS_IGNORE);
463#else
464 // Do nothing.
465#endif
466#endif
467 }
468 bool r_val = false;
469 if (sweep == 1) {
470 work_[2] = work_[0];
471 } else {
472 r_val = (work_[0] < tol2 * work_[2]);
473 }
474
475 // adjust sweep step
476 const auto adjusted_sweep_step = 2 * sweep_step_;
477 if (adjusted_sweep_step < sweep_step_upper_bound_) {
478 sweep_step_ = adjusted_sweep_step;
479 } else {
480 sweep_step_ = sweep_step_upper_bound_;
481 }
482 return r_val;
483 }
484
485 // After termination has occurred, finalize the norms for use in
486 // get_norms{0,final}.
487 void finalize() {
488 work_[0] = std::sqrt(work_[0]); // after converged
489 if (work_[2] >= 0)
490 work_[2] = std::sqrt(work_[2]); // first norm
491 // if work_[2] is minus one, then norm is not requested.
492 }
493
494 // Report norms to the caller.
495 const magnitude_type getNorms0() const { return work_[2]; }
496 const magnitude_type getNormsFinal() const { return work_[0]; }
497};
498
499template <typename MatrixType>
500void reduceVector(const ConstUnmanaged<typename BlockHelperDetails::ImplType<MatrixType>::impl_scalar_type_1d_view> zz,
501 /* */ typename BlockHelperDetails::ImplType<MatrixType>::magnitude_type *vals) {
502 IFPACK2_BLOCKHELPER_PROFILER_REGION_BEGIN;
503 IFPACK2_BLOCKHELPER_TIMER("BlockTriDi::ReduceVector", ReduceVector);
504
505 using impl_type = BlockHelperDetails::ImplType<MatrixType>;
506 using local_ordinal_type = typename impl_type::local_ordinal_type;
507 using impl_scalar_type = typename impl_type::impl_scalar_type;
508#if 0
509 const auto norm2 = KokkosBlas::nrm1(zz);
510#else
511 impl_scalar_type norm2(0);
512 Kokkos::parallel_reduce(
513 "ReduceMultiVector::Device",
514 Kokkos::RangePolicy<typename impl_type::execution_space>(0, zz.extent(0)),
515 KOKKOS_LAMBDA(const local_ordinal_type &i, impl_scalar_type &update) {
516 update += zz(i);
517 },
518 norm2);
519#endif
520 vals[0] = KokkosKernels::ArithTraits<impl_scalar_type>::abs(norm2);
521
522 IFPACK2_BLOCKHELPER_PROFILER_REGION_END;
523 IFPACK2_BLOCKHELPER_TIMER_FENCE(typename ImplType<MatrixType>::execution_space)
524}
525
526} // namespace BlockHelperDetails
527
528} // namespace Ifpack2
529
530#endif
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Definition Ifpack2_BlockHelper.hpp:211
Definition Ifpack2_BlockHelper.hpp:84
Definition Ifpack2_BlockHelper.hpp:142
Definition Ifpack2_BlockHelper.hpp:270
KB::Vector< T, l > Vector
Definition Ifpack2_BlockHelper.hpp:329
size_t size_type
Definition Ifpack2_BlockHelper.hpp:274
node_type::device_type node_device_type
Definition Ifpack2_BlockHelper.hpp:297
Kokkos::DefaultHostExecutionSpace host_execution_space
Definition Ifpack2_BlockHelper.hpp:292
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
Definition Ifpack2_BlockHelper.hpp:377
Definition Ifpack2_BlockHelper.hpp:236
Definition Ifpack2_BlockHelper.hpp:68
Definition Ifpack2_BlockHelper.hpp:95
Definition Ifpack2_BlockHelper.hpp:109
Definition Ifpack2_BlockHelper.hpp:123