10 #ifndef IFPACK2_DETAILS_CHEBYSHEV_DECL_HPP
11 #define IFPACK2_DETAILS_CHEBYSHEV_DECL_HPP
20 #include "Ifpack2_ConfigDefs.hpp"
22 #include "Teuchos_Describable.hpp"
23 #include "Tpetra_CrsMatrix.hpp"
28 #ifndef DOXYGEN_SHOULD_SKIP_THIS
29 template<
class TpetraOperatorType>
30 class ChebyshevKernel;
31 #endif // DOXYGEN_SHOULD_SKIP_THIS
74 template<
class ScalarType,
class MV>
81 typedef ScalarType
ST;
87 typedef Tpetra::Operator<
typename MV::scalar_type,
88 typename MV::local_ordinal_type,
89 typename MV::global_ordinal_type,
92 typedef Tpetra::RowMatrix<
typename MV::scalar_type,
93 typename MV::local_ordinal_type,
94 typename MV::global_ordinal_type,
97 typedef Tpetra::Vector<
typename MV::scalar_type,
98 typename MV::local_ordinal_type,
99 typename MV::global_ordinal_type,
100 typename MV::node_type>
V;
102 typedef Tpetra::Map<
typename MV::local_ordinal_type,
103 typename MV::global_ordinal_type,
237 void setZeroStartingSolution (
bool zeroStartingSolution) { zeroStartingSolution_ = zeroStartingSolution; }
292 ST getLambdaMaxForApply()
const;
308 void print (std::ostream& out);
351 typedef Kokkos::View<size_t*, typename MV::node_type::device_type> offsets_type;
358 offsets_type diagOffsets_;
367 bool savedDiagOffsets_;
384 ST computedLambdaMax_;
391 ST computedLambdaMin_;
399 ST lambdaMaxForApply_;
410 ST lambdaMinForApply_;
413 ST eigRatioForApply_;
456 bool eigKeepVectors_;
459 std::string eigenAnalysisType_;
466 int eigNormalizationFreq_;
469 bool zeroStartingSolution_;
477 bool assumeMatrixUnchanged_;
480 std::string chebyshevAlgorithm_;
483 bool computeMaxResNorm_;
486 bool computeSpectralRadius_;
490 bool ckUseNativeSpMV_;
505 void checkConstructorInput ()
const;
508 void checkInputMatrix ()
const;
536 firstIterationWithZeroStartingSolution
538 const ScalarType& alpha,
545 computeResidual (MV& R,
const MV& B,
const op_type& A,
const MV& X,
553 static void solve (MV& Z,
const V& D_inv,
const MV& R);
560 static void solve (MV& Z,
const ST alpha,
const V& D_inv,
const MV& R);
570 makeInverseDiagonal (
const row_matrix_type& A,
const bool useDiagOffsets=
false)
const;
606 textbookApplyImpl (
const op_type& A,
613 const V& D_inv)
const;
637 fourthKindApplyImpl (
const op_type& A,
667 ifpackApplyImpl (
const op_type& A,
689 cgMethodWithInitGuess (
const op_type& A,
const V& D_inv,
const int numIters,
V& x);
701 cgMethod (
const op_type& A,
const V& D_inv,
const int numIters);
704 static MT maxNormInf (
const MV& X);
742 const ST zero = Teuchos::as<ST> (0);
743 const ST one = Teuchos::as<ST> (1);
744 const ST two = Teuchos::as<ST> (2);
746 MV pAux (B.getMap (), B.getNumVectors ());
747 MV dk (B.getMap (), B.getNumVectors ());
748 MV R (B.getMap (), B.getNumVectors ());
750 ST beta = Teuchos::as<ST> (1.1) * lambdaMax;
751 ST alpha = lambdaMax / eigRatio;
753 ST delta = (beta - alpha) / two;
754 ST theta = (beta + alpha) / two;
755 ST s1 = theta / delta;
766 if (lambdaMin == one && lambdaMin == lambdaMax) {
775 if (! zeroStartingSolution_) {
779 R.update (-one, pAux, one);
780 dk.elementWiseMultiply (one/theta, D_inv, R, zero);
781 X.update (one, dk, one);
783 dk.elementWiseMultiply (one/theta, D_inv, B, zero);
787 ST rhokp1, dtemp1, dtemp2;
788 for (
int k = 0; k < numIters-1; ++k) {
790 rhokp1 = one / (two*s1 - rhok);
791 dtemp1 = rhokp1*rhok;
792 dtemp2 = two*rhokp1/delta;
796 R.update (-one, pAux, one);
798 dk.elementWiseMultiply (dtemp2, D_inv, B, dtemp1);
799 X.update (one, dk, one);
809 #endif // IFPACK2_DETAILS_CHEBYSHEV_DECL_HPP
void setParameters(Teuchos::ParameterList &plist)
Set (or reset) parameters.
Definition: Ifpack2_Details_Chebyshev_def.hpp:349
MT apply(const MV &B, MV &X)
Solve Ax=b for x with Chebyshev iteration with left diagonal scaling.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1010
Tpetra::Operator< typename MV::scalar_type, typename MV::local_ordinal_type, typename MV::global_ordinal_type, typename MV::node_type > op_type
Specialization of Tpetra::Operator.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:90
scalar_type magnitudeType
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print a description of the Chebyshev solver to out.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1703
void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Set the matrix.
Definition: Ifpack2_Details_Chebyshev_def.hpp:776
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the matrix given to the constructor.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1634
Chebyshev(Teuchos::RCP< const row_matrix_type > A)
Definition: Ifpack2_Details_Chebyshev_def.hpp:283
Tpetra::Map< typename MV::local_ordinal_type, typename MV::global_ordinal_type, typename MV::node_type > map_type
Specialization of Tpetra::Map.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:104
void print(std::ostream &out)
Print instance data to the given output stream.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1066
STS::magnitudeType MT
The type of the absolute value of a ScalarType.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:85
ScalarType ST
The type of entries in the matrix and vectors.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:81
Tpetra::RowMatrix< typename MV::scalar_type, typename MV::local_ordinal_type, typename MV::global_ordinal_type, typename MV::node_type > row_matrix_type
Specialization of Tpetra::RowMatrix.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:95
Left-scaled Chebyshev iteration.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:75
static const EVerbosityLevel verbLevel_default
std::string description() const
A single-line description of the Chebyshev solver.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1682
Teuchos::ScalarTraits< ScalarType > STS
Traits class for ST.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:83
bool hasTransposeApply() const
Whether it's possible to apply the transpose of this operator.
Definition: Ifpack2_Details_Chebyshev_def.hpp:1641
void compute()
(Re)compute the left scaling D_inv, and estimate min and max eigenvalues of D_inv * A...
Definition: Ifpack2_Details_Chebyshev_def.hpp:817
Tpetra::Vector< typename MV::scalar_type, typename MV::local_ordinal_type, typename MV::global_ordinal_type, typename MV::node_type > V
Specialization of Tpetra::Vector.
Definition: Ifpack2_Details_Chebyshev_decl.hpp:100