Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_BlockTriDiContainer_decl.hpp
Go to the documentation of this file.
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_BLOCKTRIDICONTAINER_DECL_HPP
11#define IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
12
15
16#include "Ifpack2_config.h"
17#include "Ifpack2_Container.hpp"
18#include "Tpetra_MultiVector.hpp"
19#include "Tpetra_Map.hpp"
20#include "Tpetra_RowMatrix.hpp"
21#include "Tpetra_BlockCrsMatrix_decl.hpp"
22#include <type_traits>
23#include <string>
24
25namespace Ifpack2 {
26
74
78namespace BlockTriDiContainerDetails {
83struct ImplSimdTag {};
84struct ImplSacadoTag {};
85
86template <typename T>
87struct ImplTag { typedef ImplNotAvailTag type; };
88template <>
89struct ImplTag<float> { typedef ImplSimdTag type; };
90template <>
91struct ImplTag<double> { typedef ImplSimdTag type; };
92template <>
93struct ImplTag<std::complex<float> > { typedef ImplSimdTag type; };
94template <>
95struct ImplTag<std::complex<double> > { typedef ImplSimdTag type; };
96
98template <typename MatrixType>
99struct ImplObject;
100} // namespace BlockTriDiContainerDetails
101
105template <typename MatrixType,
106 typename ImplTagType = typename BlockTriDiContainerDetails::ImplTag<typename MatrixType::scalar_type>::type>
108
114template <typename MatrixType>
115class BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
116 : public Container<MatrixType> {
118
119 private:
125 typedef MatrixType matrix_type;
126
128 typedef typename MatrixType::scalar_type scalar_type;
130 typedef typename KokkosKernels::ArithTraits<scalar_type>::magnitudeType magnitude_type;
132 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
134 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
136 typedef typename Container<MatrixType>::node_type node_type;
137
138 typedef typename Container<MatrixType>::mv_type mv_type;
139 typedef typename Container<MatrixType>::map_type map_type;
140 typedef typename Container<MatrixType>::vector_type vector_type;
141 typedef typename Container<MatrixType>::import_type import_type;
142
143 typedef typename Container<MatrixType>::HostView host_view_type;
144 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
145 typedef host_view_type HostView;
146 typedef const_host_view_type ConstHostView;
147 // typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
148 // typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
149
150 typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> crs_matrix_type;
151 typedef Tpetra::BlockCrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> block_crs_matrix_type;
152
153 const Teuchos::Array<Teuchos::Array<local_ordinal_type> > partitions_;
154
161 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
162
163 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
164 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
166 public:
168
169
185 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
186 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
187 const Teuchos::RCP<const import_type>& importer,
188 bool pointIndexed);
189
204 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
205 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
206 const int n_subparts_per_part = 1,
207 bool overlapCommAndComp = false,
208 bool useSequentialMethod = false,
209 const int block_size = -1,
210 const bool explicitConversion = false);
211
213 ~BlockTriDiContainer() override;
214
215 struct ComputeParameters {
222 magnitude_type addRadiallyToDiagonal = KokkosKernels::ArithTraits<magnitude_type>::zero();
223 };
224
226 struct ApplyParameters {
229 bool zeroStartingSolution = false;
231 scalar_type dampingFactor = KokkosKernels::ArithTraits<scalar_type>::one();
234 int maxNumSweeps = 1;
244 magnitude_type tolerance = KokkosKernels::ArithTraits<magnitude_type>::zero();
252 int checkToleranceEvery = 1;
253 };
254
256
258
260 void setParameters(const Teuchos::ParameterList& List) override;
261
262 void clearBlocks() override;
263
265
267
269 void initialize() override;
270
272 void compute() override;
273
274 // \brief Compute <tt>Y := D^{-1} (X - R*Y)</tt>.
275 void applyInverseJacobi(const mv_type& X, mv_type& Y,
276 scalar_type dampingFactor,
277 bool zeroStartingSolution = false,
278 int numSweeps = 1) const override;
279
281 ComputeParameters createDefaultComputeParameters() const;
282
294 void compute(const ComputeParameters& input);
295
297 ApplyParameters createDefaultApplyParameters() const;
298
305 int applyInverseJacobi(const mv_type& X, mv_type& Y,
306 const ApplyParameters& input) const;
307
311 const magnitude_type getNorms0() const;
312
315 const magnitude_type getNormsFinal() const;
316
319 void
320 apply(const_host_view_type X,
321 host_view_type Y,
322 int blockIndex,
323 Teuchos::ETransp mode = Teuchos::NO_TRANS,
324 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
325 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
326
329 void
330 weightedApply(const_host_view_type X,
331 host_view_type Y,
332 const_host_view_type W,
333 int blockIndex,
334 Teuchos::ETransp mode = Teuchos::NO_TRANS,
335 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
336 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
337
339
341
345 std::ostream& print(std::ostream& os) const override;
346
348
350
352 std::string description() const override;
353
355 void
356 describe(Teuchos::FancyOStream& out,
357 const Teuchos::EVerbosityLevel verbLevel =
358 Teuchos::Describable::verbLevel_default) const override;
359
361
363 static std::string getName();
364
365 private:
368
369 // hide details of impl using ImplObj; finally I understand why AMB did that way.
370 Teuchos::RCP<BlockTriDiContainerDetails::ImplObject<MatrixType> > impl_;
371 int n_subparts_per_part_;
372 int block_size_ = -1;
373
374 // initialize distributed and local objects
375 void initInternal(const Teuchos::RCP<const row_matrix_type>& matrix,
376 const Teuchos::RCP<const import_type>& importer,
377 const bool overlapCommAndComp,
378 const bool useSeqMethod,
379 const int block_size = -1,
380 const bool explicitConversion = false);
381
382 void clearInternal();
383
384 // Decide whether the fused block Jacobi path can and should be used.
385 bool shouldUseFusedBlockJacobi(
386 const Teuchos::RCP<const row_matrix_type>& matrix,
387 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
388 bool useSeqMethod);
389};
390
398template <typename MatrixType>
399class BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag>
400 : public Container<MatrixType> {
401 private:
402 typedef typename MatrixType::scalar_type scalar_type;
403 typedef typename KokkosKernels::ArithTraits<scalar_type>::magnitudeType magnitude_type;
404 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
405 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
406
407 typedef typename Container<MatrixType>::mv_type mv_type;
408 typedef typename Container<MatrixType>::import_type import_type;
409
410 typedef typename Container<MatrixType>::HostView host_view_type;
411 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
412 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
413
414 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
415 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
416
417 public:
418 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
419 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
420 const Teuchos::RCP<const import_type>& importer,
421 bool pointIndexed)
422 : Container<MatrixType>(matrix, partitions, pointIndexed) {
423 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Error: BlockTriDiContainer is not available for this scalar_type");
424 }
425
426 void setParameters(const Teuchos::ParameterList& List) override {}
427 void clearBlocks() override {}
428
429 void initialize() override {}
430 void compute() override {}
431 void applyInverseJacobi(const mv_type& X, mv_type& Y,
432 scalar_type dampingFactor,
433 bool zeroStartingSolution = false,
434 int numSweeps = 1) const override {}
435
436 void
437 apply(const_host_view_type X,
438 host_view_type Y,
439 int blockIndex,
440 Teuchos::ETransp mode = Teuchos::NO_TRANS,
441 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
442 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
443
444 void
445 weightedApply(const_host_view_type X,
446 host_view_type Y,
447 const_host_view_type W,
448 int blockIndex,
449 Teuchos::ETransp mode = Teuchos::NO_TRANS,
450 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
451 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
452
453 std::ostream& print(std::ostream& os) const override {
454 return os << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
455 }
456
457 std::string description() const override {
458 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
459 }
460
461 void
462 describe(Teuchos::FancyOStream& out,
463 const Teuchos::EVerbosityLevel verbLevel =
464 Teuchos::Describable::verbLevel_default) const override {
465 out << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
466 }
467
468 static std::string getName() {
469 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
470 }
471};
472
473} // namespace Ifpack2
474
475#endif // IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
void setParameters(const Teuchos::ParameterList &List) override
Set parameters, if any.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:426
void weightedApply(const_host_view_type X, host_view_type Y, const_host_view_type W, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:445
void compute() override
Extract the local diagonal blocks and prepare the solver.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:430
void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:429
void apply(const_host_view_type X, host_view_type Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
Compute Y := alpha * M^{-1} X + beta*Y.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:437
void applyInverseJacobi(const mv_type &X, mv_type &Y, scalar_type dampingFactor, bool zeroStartingSolution=false, int numSweeps=1) const override
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).
Definition Ifpack2_BlockTriDiContainer_decl.hpp:431
std::ostream & print(std::ostream &os) const override
Print basic information about the container to os.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:453
Store and solve local block tridiagonal linear problems.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:107
Interface for creating and solving a set of local linear problems.
Definition Ifpack2_Container_decl.hpp:79
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:105
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Definition Ifpack2_BlockTriDiContainer_decl.hpp:82