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 "Ifpack2_BlockHelper_ETI.hpp"
23#include <type_traits>
24#include <string>
25
26namespace Ifpack2 {
27
75
79template <typename MatrixType,
80 typename ImplTagType = typename BlockTriDiContainerDetails::ImplTag<typename MatrixType::scalar_type>::type>
82
88template <typename MatrixType>
89class BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
90 : public Container<MatrixType> {
92
93 private:
99 typedef MatrixType matrix_type;
100
102 typedef typename MatrixType::scalar_type scalar_type;
104 typedef typename KokkosKernels::ArithTraits<scalar_type>::magnitudeType magnitude_type;
106 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
108 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
110 typedef typename Container<MatrixType>::node_type node_type;
111
112 typedef typename Container<MatrixType>::mv_type mv_type;
113 typedef typename Container<MatrixType>::map_type map_type;
114 typedef typename Container<MatrixType>::vector_type vector_type;
115 typedef typename Container<MatrixType>::import_type import_type;
116
117 typedef typename Container<MatrixType>::HostView host_view_type;
118 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
119 typedef host_view_type HostView;
120 typedef const_host_view_type ConstHostView;
121 // typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
122 // typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
123
124 typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> crs_matrix_type;
125 typedef Tpetra::BlockCrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> block_crs_matrix_type;
126
127 const Teuchos::Array<Teuchos::Array<local_ordinal_type> > partitions_;
128
135 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
136
137 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
138 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
140 public:
142
143
159 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
160 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
161 const Teuchos::RCP<const import_type>& importer,
162 bool pointIndexed);
163
178 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
179 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
180 const int n_subparts_per_part = 1,
181 bool overlapCommAndComp = false,
182 bool useSequentialMethod = false,
183 const int block_size = -1,
184 const bool explicitConversion = false);
185
187 ~BlockTriDiContainer() override;
188
189 struct ComputeParameters {
196 magnitude_type addRadiallyToDiagonal = KokkosKernels::ArithTraits<magnitude_type>::zero();
197 };
198
200 struct ApplyParameters {
203 bool zeroStartingSolution = false;
205 scalar_type dampingFactor = KokkosKernels::ArithTraits<scalar_type>::one();
208 int maxNumSweeps = 1;
218 magnitude_type tolerance = KokkosKernels::ArithTraits<magnitude_type>::zero();
226 int checkToleranceEvery = 1;
227 };
228
230
232
234 void setParameters(const Teuchos::ParameterList& List) override;
235
236 void clearBlocks() override;
237
239
241
243 void initialize() override;
244
246 void compute() override;
247
248 // \brief Compute <tt>Y := D^{-1} (X - R*Y)</tt>.
249 void applyInverseJacobi(const mv_type& X, mv_type& Y,
250 scalar_type dampingFactor,
251 bool zeroStartingSolution = false,
252 int numSweeps = 1) const override;
253
255 ComputeParameters createDefaultComputeParameters() const;
256
268 void compute(const ComputeParameters& input);
269
271 ApplyParameters createDefaultApplyParameters() const;
272
279 int applyInverseJacobi(const mv_type& X, mv_type& Y,
280 const ApplyParameters& input) const;
281
285 const magnitude_type getNorms0() const;
286
289 const magnitude_type getNormsFinal() const;
290
293 void
294 apply(const_host_view_type X,
295 host_view_type Y,
296 int blockIndex,
297 Teuchos::ETransp mode = Teuchos::NO_TRANS,
298 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
299 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
300
303 void
304 weightedApply(const_host_view_type X,
305 host_view_type Y,
306 const_host_view_type W,
307 int blockIndex,
308 Teuchos::ETransp mode = Teuchos::NO_TRANS,
309 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
310 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
311
313
315
319 std::ostream& print(std::ostream& os) const override;
320
322
324
326 std::string description() const override;
327
329 void
330 describe(Teuchos::FancyOStream& out,
331 const Teuchos::EVerbosityLevel verbLevel =
332 Teuchos::Describable::verbLevel_default) const override;
333
335
337 static std::string getName();
338
339 private:
342
343 // hide details of impl using ImplObj; finally I understand why AMB did that way.
344 Teuchos::RCP<BlockTriDiContainerDetails::ImplObject<MatrixType> > impl_;
345 int n_subparts_per_part_;
346 int block_size_ = -1;
347
348 // initialize distributed and local objects
349 void initInternal(const Teuchos::RCP<const row_matrix_type>& matrix,
350 const Teuchos::RCP<const import_type>& importer,
351 const bool overlapCommAndComp,
352 const bool useSeqMethod,
353 const int block_size = -1,
354 const bool explicitConversion = false);
355
356 void clearInternal();
357
358 // Decide whether the fused block Jacobi path can and should be used.
359 bool shouldUseFusedBlockJacobi(
360 const Teuchos::RCP<const row_matrix_type>& matrix,
361 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
362 bool useSeqMethod);
363};
364
372template <typename MatrixType>
373class BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag>
374 : public Container<MatrixType> {
375 private:
376 typedef typename MatrixType::scalar_type scalar_type;
377 typedef typename KokkosKernels::ArithTraits<scalar_type>::magnitudeType magnitude_type;
378 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
379 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
380
381 typedef typename Container<MatrixType>::mv_type mv_type;
382 typedef typename Container<MatrixType>::import_type import_type;
383
384 typedef typename Container<MatrixType>::HostView host_view_type;
385 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
386 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
387
388 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
389 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
390
391 public:
392 BlockTriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
393 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
394 const Teuchos::RCP<const import_type>& importer,
395 bool pointIndexed)
396 : Container<MatrixType>(matrix, partitions, pointIndexed) {
397 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Error: BlockTriDiContainer is not available for this scalar_type");
398 }
399
400 void setParameters(const Teuchos::ParameterList& List) override {}
401 void clearBlocks() override {}
402
403 void initialize() override {}
404 void compute() override {}
405 void applyInverseJacobi(const mv_type& X, mv_type& Y,
406 scalar_type dampingFactor,
407 bool zeroStartingSolution = false,
408 int numSweeps = 1) const override {}
409
410 void
411 apply(const_host_view_type X,
412 host_view_type Y,
413 int blockIndex,
414 Teuchos::ETransp mode = Teuchos::NO_TRANS,
415 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
416 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
417
418 void
419 weightedApply(const_host_view_type X,
420 host_view_type Y,
421 const_host_view_type W,
422 int blockIndex,
423 Teuchos::ETransp mode = Teuchos::NO_TRANS,
424 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
425 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
426
427 std::ostream& print(std::ostream& os) const override {
428 return os << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
429 }
430
431 std::string description() const override {
432 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
433 }
434
435 void
436 describe(Teuchos::FancyOStream& out,
437 const Teuchos::EVerbosityLevel verbLevel =
438 Teuchos::Describable::verbLevel_default) const override {
439 out << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
440 }
441
442 static std::string getName() {
443 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
444 }
445};
446
447} // namespace Ifpack2
448
449#endif // IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
void setParameters(const Teuchos::ParameterList &List) override
Set parameters, if any.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:400
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:419
void compute() override
Extract the local diagonal blocks and prepare the solver.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:404
void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:403
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:411
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:405
std::ostream & print(std::ostream &os) const override
Print basic information about the container to os.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:427
Store and solve local block tridiagonal linear problems.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:81
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