Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_BlockRelaxation_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_BLOCKRELAXATION_DECL_HPP
11#define IFPACK2_BLOCKRELAXATION_DECL_HPP
12
15
17#include "Ifpack2_Partitioner.hpp"
19#include "Ifpack2_ContainerFactory.hpp"
20#include "Tpetra_BlockCrsMatrix.hpp"
21#include <type_traits>
22
23namespace Ifpack2 {
24
48template <class MatrixType, class ContainerType = Container<MatrixType>>
49class BlockRelaxation : virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
50 typename MatrixType::local_ordinal_type,
51 typename MatrixType::global_ordinal_type,
52 typename MatrixType::node_type>,
53 virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
54 typename MatrixType::local_ordinal_type,
55 typename MatrixType::global_ordinal_type,
56 typename MatrixType::node_type>> {
57 public:
59
60
62 typedef typename MatrixType::scalar_type scalar_type;
63
65 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
66
68 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
70 typedef typename MatrixType::node_type node_type;
71
73 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
74
76 typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> row_matrix_type;
77
78 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
79 "Ifpack2::BlockRelaxation: Please use MatrixType = Tpetra::RowMatrix.");
80
81 static_assert(std::is_same<ContainerType, Container<row_matrix_type>>::value,
82 "Ifpack2::BlockRelaxation: Do NOT specify the (second) "
83 "ContainerType template parameter explicitly. The default "
84 "value is fine. Please instead specify the container type to "
85 "use by setting the \"relaxation: container\" parameter.");
86
88 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
89
90 private:
95 void setParametersImpl(Teuchos::ParameterList& params);
96
97 void computeImporter() const;
98
100
101 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type,
103 MV;
104 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vector_type;
105 typedef Teuchos::ScalarTraits<scalar_type> STS;
106 typedef Teuchos::ScalarTraits<magnitude_type> STM;
107 typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type,
109 crs_matrix_type;
110 typedef Tpetra::BlockCrsMatrix<scalar_type, local_ordinal_type,
112 block_crs_matrix_type;
113 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
114
115 public:
117
118 // \name Constructors and Destructors
120
150 explicit BlockRelaxation(const Teuchos::RCP<const row_matrix_type>& Matrix);
151
153 virtual ~BlockRelaxation();
154
156
158
160
187 void setParameters(const Teuchos::ParameterList& params);
188
189 bool supportsZeroStartingSolution() { return true; }
190
191 void setZeroStartingSolution(bool zeroStartingSolution) { ZeroStartingSolution_ = zeroStartingSolution; };
192
194 Teuchos::RCP<const Teuchos::ParameterList>
195 getValidParameters() const;
196
198 void initialize();
199
201 inline bool isInitialized() const {
202 return (IsInitialized_);
203 }
204
206 void compute();
207
209 inline bool isComputed() const {
210 return (IsComputed_);
211 }
212
214
216
239 virtual void
240 setMatrix(const Teuchos::RCP<const row_matrix_type>& A);
241
243
245
247
257 void apply(const MV& X,
258 MV& Y,
259 Teuchos::ETransp mode = Teuchos::NO_TRANS,
260 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
261 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
262
264 Teuchos::RCP<const map_type> getDomainMap() const;
265
267 Teuchos::RCP<const map_type> getRangeMap() const;
268
269 bool hasTransposeApply() const;
270
272
278 void applyMat(const MV& X,
279 MV& Y,
280 Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
281
283
285
287 Teuchos::RCP<const Teuchos::Comm<int>> getComm() const;
288
290 Teuchos::RCP<const row_matrix_type> getMatrix() const;
291
293 int getNumInitialize() const;
294
296 int getNumCompute() const;
297
299 int getNumApply() const;
300
302 double getInitializeTime() const;
303
305 double getComputeTime() const;
306
308 double getApplyTime() const;
309
311 size_t getNodeSmootherComplexity() const;
312
314
316
318 std::string description() const;
319
321 void
322 describe(Teuchos::FancyOStream& out,
323 const Teuchos::EVerbosityLevel verbLevel =
324 Teuchos::Describable::verbLevel_default) const;
325
327
329 Teuchos::RCP<Ifpack2::Partitioner<Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>>> getPartitioner() { return Partitioner_; }
330
331 private:
334
337 operator=(const BlockRelaxation<MatrixType, ContainerType>& RHS);
338
339 virtual void ApplyInverseJacobi(const MV& X, MV& Y) const;
340
341 virtual void ApplyInverseGS(const MV& X, MV& Y) const;
342
343 virtual void ApplyInverseSGS(const MV& X, MV& Y) const;
344
347 void ExtractSubmatricesStructure();
348
350
352
354 Teuchos::RCP<const row_matrix_type> A_;
355
357 mutable Teuchos::RCP<Container<row_matrix_type>> Container_;
358
359 // FIXME (mfh 06 Oct 2014) This doesn't comply with the naming
360 // convention for instance members of a class. Furthermore, the
361 // class should keep the Vector, not the ArrayRCP to the data _in_
362 // the Vector.
363 // FIXED! (amk 10 Nov 2015)
364 mutable Teuchos::RCP<vector_type> DiagRCP_;
365
367 Teuchos::RCP<Ifpack2::Partitioner<Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>>> Partitioner_;
368
371 std::string PartitionerType_;
372
374 Teuchos::ParameterList List_;
375
377 int NumSweeps_;
378
380 local_ordinal_type NumLocalBlocks_;
381
383 std::string containerType_;
384
386 Details::RelaxationType PrecType_;
387
389 bool IsParallel_;
390
392 bool ZeroStartingSolution_;
393
396 bool hasBlockCrsMatrix_;
397
399 bool DoBackwardGS_;
400
402 int OverlapLevel_;
403
405 // only valid with block Jacobi relaxation and overlapping blocks (e.g., defined by "partitioner: type" "user" and "parts: " or "global ID parts:"). Average solutions in overlapped regions (i.e., after summing different solutions divide by number of blocks contain this dof). When false (the default) symmetric averaging performed (i.e., average residuals and solutions).
406 bool nonsymCombine_;
407
409 std::string schwarzCombineMode_;
410
412 scalar_type DampingFactor_;
413
415 bool decouple_;
416
418 bool IsInitialized_;
419
421 bool IsComputed_;
422
424 int NumInitialize_;
425
427 int NumCompute_;
428
430 bool TimerForApply_;
431
433 mutable int NumApply_;
434
436 double InitializeTime_;
437
439 double ComputeTime_;
440
442 mutable double ApplyTime_;
443
445 local_ordinal_type NumLocalRows_;
446
448 global_ordinal_type NumGlobalRows_;
449
451 global_ordinal_type NumGlobalNonzeros_;
452
455 Teuchos::RCP<vector_type> W_;
456
457 mutable Teuchos::RCP<const Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type>> Importer_;
458
460}; // class BlockRelaxation
461
462} // namespace Ifpack2
463
464#endif // IFPACK2_BLOCKRELAXATION_DECL_HPP
Declaration of interface for preconditioners that can change their matrix after construction.
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition Ifpack2_BlockRelaxation_decl.hpp:56
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner's constructor.
Definition Ifpack2_BlockRelaxation_def.hpp:340
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization corresponding to MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:76
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:349
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner's parameters.
Definition Ifpack2_BlockRelaxation_decl.hpp:191
int getNumApply() const
Returns the number of calls to apply().
Definition Ifpack2_BlockRelaxation_def.hpp:390
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition Ifpack2_BlockRelaxation_def.hpp:96
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition Ifpack2_BlockRelaxation_def.hpp:742
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_BlockRelaxation_def.hpp:416
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
Definition Ifpack2_BlockRelaxation_def.hpp:1068
Tpetra::Import< local_ordinal_type, global_ordinal_type, node_type > import_type
Tpetra::Importer specialization for use with MatrixType and compatible MultiVectors.
Definition Ifpack2_BlockRelaxation_decl.hpp:88
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:65
virtual ~BlockRelaxation()
Destructor.
Definition Ifpack2_BlockRelaxation_def.hpp:91
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition Ifpack2_BlockRelaxation_decl.hpp:201
int getNumCompute() const
Returns the number of calls to compute().
Definition Ifpack2_BlockRelaxation_def.hpp:384
double getApplyTime() const
Returns the time spent in apply().
Definition Ifpack2_BlockRelaxation_def.hpp:411
int getNumInitialize() const
Returns the number of calls to initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:378
MatrixType::node_type node_type
Node type of the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:70
void initialize()
Initialize.
Definition Ifpack2_BlockRelaxation_def.hpp:555
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:362
bool isComputed() const
Return true if compute() has been called.
Definition Ifpack2_BlockRelaxation_decl.hpp:209
void applyMat(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS) const
Applies the matrix to a Tpetra::MultiVector.
Definition Ifpack2_BlockRelaxation_def.hpp:541
double getComputeTime() const
Returns the time spent in compute().
Definition Ifpack2_BlockRelaxation_def.hpp:404
void setParameters(const Teuchos::ParameterList &params)
Sets all the parameters for the preconditioner.
Definition Ifpack2_BlockRelaxation_def.hpp:156
double getInitializeTime() const
Returns the time spent in initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:397
void apply(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Applies the preconditioner to X, returns the result in Y.
Definition Ifpack2_BlockRelaxation_def.hpp:433
std::string description() const
A one-line description of this object.
Definition Ifpack2_BlockRelaxation_def.hpp:1007
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:68
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition Ifpack2_BlockRelaxation_def.hpp:327
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_BlockRelaxation_decl.hpp:73
Teuchos::RCP< Ifpack2::Partitioner< Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > > getPartitioner()
For diagnostic purposes.
Definition Ifpack2_BlockRelaxation_decl.hpp:329
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_BlockRelaxation_def.hpp:31
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:62
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40