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
78 namespace BlockTriDiContainerDetails {
82 struct ImplNotAvailTag {};
83 struct ImplSimdTag {};
84 struct ImplSacadoTag {};
85
86 template<typename T> struct ImplTag { typedef ImplNotAvailTag type; };
87 template<> struct ImplTag<float> { typedef ImplSimdTag type; };
88 template<> struct ImplTag<double> { typedef ImplSimdTag type; };
89 template<> struct ImplTag<std::complex<float> > { typedef ImplSimdTag type; };
90 template<> struct ImplTag<std::complex<double> > { typedef ImplSimdTag type; };
91
93 template<typename MatrixType> struct ImplObject;
94 }
95
99 template <typename MatrixType,
100 typename ImplTagType = typename BlockTriDiContainerDetails::ImplTag<typename MatrixType::scalar_type>::type>
102
108 template <typename MatrixType>
109 class BlockTriDiContainer<MatrixType,BlockTriDiContainerDetails::ImplSimdTag>
110 : public Container<MatrixType> {
112
113 private:
119 typedef MatrixType matrix_type;
120
122 typedef typename MatrixType::scalar_type scalar_type;
124 typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
126 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
128 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
130 typedef typename Container<MatrixType>::node_type node_type;
131
132 typedef typename Container<MatrixType>::mv_type mv_type;
133 typedef typename Container<MatrixType>::map_type map_type;
134 typedef typename Container<MatrixType>::vector_type vector_type;
135 typedef typename Container<MatrixType>::import_type import_type;
136
137 typedef typename Container<MatrixType>::HostView host_view_type;
138 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
139 typedef host_view_type HostView;
140 typedef const_host_view_type ConstHostView;
141 //typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
142 //typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
143
144 typedef Tpetra::CrsMatrix
145 <scalar_type,local_ordinal_type,global_ordinal_type,node_type> crs_matrix_type;
146 typedef Tpetra::BlockCrsMatrix
147 <scalar_type,local_ordinal_type,global_ordinal_type,node_type> block_crs_matrix_type;
148
149 const Teuchos::Array<Teuchos::Array<local_ordinal_type> > partitions_;
150
157 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
158
159 static_assert (std::is_same<MatrixType, row_matrix_type>::value,
160 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
162 public:
164
165
181 BlockTriDiContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
182 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
183 const Teuchos::RCP<const import_type>& importer,
184 bool pointIndexed);
185
200 BlockTriDiContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
201 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
202 const int n_subparts_per_part = 1,
203 bool overlapCommAndComp = false,
204 bool useSequentialMethod = false,
205 const int block_size = -1,
206 const bool explicitConversion = false);
207
209 ~BlockTriDiContainer () override;
210
211 struct ComputeParameters {
218 magnitude_type addRadiallyToDiagonal = Kokkos::ArithTraits<magnitude_type>::zero();
219 };
220
222 struct ApplyParameters {
225 bool zeroStartingSolution = false;
227 scalar_type dampingFactor = Kokkos::ArithTraits<scalar_type>::one();
230 int maxNumSweeps = 1;
240 magnitude_type tolerance = Kokkos::ArithTraits<magnitude_type>::zero();
248 int checkToleranceEvery = 1;
249 };
250
252
254
256 void setParameters(const Teuchos::ParameterList& List) override;
257
258 void clearBlocks() override;
259
261
263
265 void initialize () override;
266
268 void compute () override;
269
270 // \brief Compute <tt>Y := D^{-1} (X - R*Y)</tt>.
271 void applyInverseJacobi (const mv_type& X, mv_type& Y,
272 scalar_type dampingFactor,
273 bool zeroStartingSolution = false,
274 int numSweeps = 1) const override;
275
277 ComputeParameters createDefaultComputeParameters () const;
278
290 void compute (const ComputeParameters& input);
291
293 ApplyParameters createDefaultApplyParameters () const;
294
301 int applyInverseJacobi (const mv_type& X, mv_type& Y,
302 const ApplyParameters& input) const;
303
307 const magnitude_type getNorms0 () const;
308
311 const magnitude_type getNormsFinal () const;
312
315 void
316 apply (const_host_view_type X,
317 host_view_type Y,
318 int blockIndex,
319 Teuchos::ETransp mode = Teuchos::NO_TRANS,
320 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
321 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
322
325 void
326 weightedApply (const_host_view_type X,
327 host_view_type Y,
328 const_host_view_type W,
329 int blockIndex,
330 Teuchos::ETransp mode = Teuchos::NO_TRANS,
331 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
332 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
333
335
337
341 std::ostream& print (std::ostream& os) const override;
342
344
346
348 std::string description () const override;
349
351 void
352 describe (Teuchos::FancyOStream &out,
353 const Teuchos::EVerbosityLevel verbLevel =
354 Teuchos::Describable::verbLevel_default) const override;
355
357
359 static std::string getName();
360
361 private:
364
365 // hide details of impl using ImplObj; finally I understand why AMB did that way.
366 Teuchos::RCP<BlockTriDiContainerDetails::ImplObject<MatrixType> > impl_;
367 int n_subparts_per_part_;
368 int block_size_ = -1;
369
370 // initialize distributed and local objects
371 void initInternal (const Teuchos::RCP<const row_matrix_type>& matrix,
372 const Teuchos::RCP<const import_type> &importer,
373 const bool overlapCommAndComp,
374 const bool useSeqMethod,
375 const int block_size = -1,
376 const bool explicitConversion = false);
377
378 void clearInternal();
379 };
380
388 template <typename MatrixType>
389 class BlockTriDiContainer<MatrixType,BlockTriDiContainerDetails::ImplNotAvailTag>
390 : public Container<MatrixType> {
391 private:
392 typedef typename MatrixType::scalar_type scalar_type;
393 typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
394 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
395 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
396
397 typedef typename Container<MatrixType>::mv_type mv_type;
398 typedef typename Container<MatrixType>::import_type import_type;
399
400 typedef typename Container<MatrixType>::HostView host_view_type;
401 typedef typename Container<MatrixType>::ConstHostView const_host_view_type;
402 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
403
404 static_assert (std::is_same<MatrixType, row_matrix_type>::value,
405 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
406 public:
407
408 BlockTriDiContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
409 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
410 const Teuchos::RCP<const import_type>& importer,
411 bool pointIndexed)
412 : Container<MatrixType>(matrix, partitions, pointIndexed) {
413 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Error: BlockTriDiContainer is not available for this scalar_type");
414 }
415
416 void setParameters(const Teuchos::ParameterList& List) override {}
417 void clearBlocks() override {}
418
419 void initialize () override {}
420 void compute () override {}
421 void applyInverseJacobi (const mv_type& X, mv_type& Y,
422 scalar_type dampingFactor,
423 bool zeroStartingSolution = false,
424 int numSweeps = 1) const override {}
425
426 void
427 apply (const_host_view_type X,
428 host_view_type Y,
429 int blockIndex,
430 Teuchos::ETransp mode = Teuchos::NO_TRANS,
431 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
432 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
433
434 void
435 weightedApply (const_host_view_type X,
436 host_view_type Y,
437 const_host_view_type W,
438 int blockIndex,
439 Teuchos::ETransp mode = Teuchos::NO_TRANS,
440 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
441 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
442
443 std::ostream& print (std::ostream& os) const override {
444 return os << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
445 }
446
447 std::string description () const override {
448 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
449 }
450
451 void
452 describe (Teuchos::FancyOStream &out,
453 const Teuchos::EVerbosityLevel verbLevel =
454 Teuchos::Describable::verbLevel_default) const override {
455 out << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
456 }
457
458 static std::string getName() {
459 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
460 }
461 };
462
463
464} // namespace Ifpack2
465
466#endif // IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
void setParameters(const Teuchos::ParameterList &List) override
Set parameters, if any.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:416
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:435
void compute() override
Extract the local diagonal blocks and prepare the solver.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:420
void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:419
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:427
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:421
std::ostream & print(std::ostream &os) const override
Print basic information about the container to os.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:443
Store and solve local block tridiagonal linear problems.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:101
Interface for creating and solving a set of local linear problems.
Definition Ifpack2_Container_decl.hpp:80
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:106
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41
Definition Ifpack2_BlockTriDiContainer_decl.hpp:82