Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Condest.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_CONDEST_HPP
11#define IFPACK2_CONDEST_HPP
12
13#include "Ifpack2_ConfigDefs.hpp"
14#include "Ifpack2_CondestType.hpp"
16#include <Teuchos_Ptr.hpp>
17
18namespace Ifpack2 {
19
59template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
60typename Teuchos::ScalarTraits<Scalar>::magnitudeType
62 const Ifpack2::CondestType CT,
63 const int MaxIters = 1550,
64 const typename Teuchos::ScalarTraits<Scalar>::magnitudeType& Tol = Teuchos::as<Scalar>(1e-9),
65 const Teuchos::Ptr<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& matrix_in = Teuchos::null) {
66 using Teuchos::Ptr;
67 typedef Teuchos::ScalarTraits<Scalar> STS;
68 typedef typename STS::magnitudeType MT;
69 typedef Teuchos::ScalarTraits<MT> STM;
70 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> row_matrix_type;
71 typedef Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> vec_type;
72
73 MT condNumEst = -STS::magnitude(STS::one());
74
75 // Users may either provide a matrix for which to estimate the
76 // condition number, or use the Preconditioner's built-in matrix.
77 Ptr<const row_matrix_type> matrix = matrix_in;
78 if (matrix_in == Teuchos::null) {
79 matrix = TIFP.getMatrix().ptr();
80 TEUCHOS_TEST_FOR_EXCEPTION(
81 matrix == Teuchos::null,
82 std::logic_error,
83 "Ifpack2::Condest: Both the input matrix (matrix_in) and the Ifpack2 "
84 "preconditioner's matrix are null, so we have no matrix with which to "
85 "compute a condition number estimate. This probably indicates a bug "
86 "in Ifpack2, since no Ifpack2::Preconditioner subclass should accept a"
87 "null matrix.");
88 }
89
90 if (CT == Ifpack2::Cheap) {
91 vec_type ones(TIFP.getDomainMap()); // Vector of ones
92 ones.putScalar(STS::one());
93 vec_type onesResult(TIFP.getRangeMap()); // A*ones
94 onesResult.putScalar(STS::zero());
95 TIFP.apply(ones, onesResult);
96 condNumEst = onesResult.normInf(); // max (abs (A*ones))
97 TEUCHOS_TEST_FOR_EXCEPTION(
98 STM::isnaninf(condNumEst),
99 std::runtime_error,
100 "Ifpack2::Condest: $\\|A*[1, ..., 1]^T\\|_{\\infty}$ = " << condNumEst << " is NaN or Inf.");
101 } else if (CT == Ifpack2::CG) {
102 TEUCHOS_TEST_FOR_EXCEPTION(
103 true, std::logic_error,
104 "Ifpack2::Condest: Condition number estimation using CG is not currently supported.");
105 } else if (CT == Ifpack2::GMRES) {
106 TEUCHOS_TEST_FOR_EXCEPTION(
107 true, std::logic_error,
108 "Ifpack2::Condest: Condition number estimation using GMRES is not currently supported.");
109 }
110 return condNumEst;
111}
112
113} // namespace Ifpack2
114
115#endif // IFPACK2_CONDEST_HPP
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
virtual Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMatrix() const =0
The input matrix given to the constructor.
virtual void apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const =0
Apply the preconditioner to X, putting the result in Y.
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const =0
The range Map of this operator.
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const =0
The domain Map of this operator.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
CondestType
Ifpack2::CondestType: enum to define the type of condition number estimate.
Definition Ifpack2_CondestType.hpp:17
@ CG
Uses AztecOO's CG.
Definition Ifpack2_CondestType.hpp:19
@ Cheap
cheap estimate
Definition Ifpack2_CondestType.hpp:18
@ GMRES
Uses AztecOO's GMRES.
Definition Ifpack2_CondestType.hpp:20