43 #ifndef IFPACK2_CONDEST_HPP 
   44 #define IFPACK2_CONDEST_HPP 
   46 #include "Ifpack2_ConfigDefs.hpp" 
   47 #include "Ifpack2_CondestType.hpp" 
   49 #include <Teuchos_Ptr.hpp> 
   92 template<
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   96          const int MaxIters = 1550,
 
   98          const Teuchos::Ptr<
const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& matrix_in = Teuchos::null)
 
  102   typedef typename STS::magnitudeType MT;
 
  104   typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> row_matrix_type;
 
  105   typedef Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> vec_type;
 
  107   MT condNumEst = -STS::magnitude( STS::one() );
 
  111   Ptr<const row_matrix_type> matrix = matrix_in;
 
  112   if (matrix_in == Teuchos::null) {
 
  115       matrix == Teuchos::null,
 
  117       "Ifpack2::Condest: Both the input matrix (matrix_in) and the Ifpack2 " 
  118       "preconditioner's matrix are null, so we have no matrix with which to " 
  119       "compute a condition number estimate.  This probably indicates a bug " 
  120       "in Ifpack2, since no Ifpack2::Preconditioner subclass should accept a" 
  126     ones.putScalar (STS::one ());
 
  128     onesResult.putScalar (STS::zero ());
 
  129     TIFP.
apply (ones, onesResult);
 
  130     condNumEst = onesResult.normInf (); 
 
  132       STM::isnaninf (condNumEst),
 
  134       "Ifpack2::Condest: $\\|A*[1, ..., 1]^T\\|_{\\infty}$ = " << condNumEst << 
" is NaN or Inf.");
 
  137       true, std::logic_error,
 
  138       "Ifpack2::Condest: Condition number estimation using CG is not currently supported.");
 
  141       true, std::logic_error,
 
  142       "Ifpack2::Condest: Condition number estimation using GMRES is not currently supported.");
 
  149 #endif // IFPACK2_CONDEST_HPP 
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const =0
The domain Map of this operator. 
 
CondestType
Ifpack2::CondestType: enum to define the type of condition number estimate. 
Definition: Ifpack2_CondestType.hpp:50
 
virtual Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMatrix() const =0
The input matrix given to the constructor. 
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const =0
The range Map of this operator. 
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
Uses AztecOO's CG. 
Definition: Ifpack2_CondestType.hpp:52
 
Uses AztecOO's GMRES. 
Definition: Ifpack2_CondestType.hpp:53
 
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. 
 
cheap estimate 
Definition: Ifpack2_CondestType.hpp:51