MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

MueLu utility class. More...

#include <MueLu_UtilitiesBase_fwd.hpp>

Inheritance diagram for MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
MueLu::Utilities< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Public Types

typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
Magnitude
 

Static Public Member Functions

static RCP< Matrix > Crs2Op (RCP< CrsMatrix > Op)
 
static RCP< CrsMatrixWrap > GetThresholdedMatrix (const RCP< Matrix > &Ain, const Scalar threshold, const bool keepDiagonal=true, const GlobalOrdinal expectedNNZperRow=-1)
 Threshold a matrix. More...
 
static RCP< Xpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > > 
GetThresholdedGraph (const RCP< Matrix > &A, const Magnitude threshold, const GlobalOrdinal expectedNNZperRow=-1)
 Threshold a graph. More...
 
static Teuchos::ArrayRCP< ScalarGetMatrixDiagonal_arcp (const Matrix &A)
 Extract Matrix Diagonal. More...
 
static RCP< Vector > GetMatrixDiagonal (const Matrix &A)
 Extract Matrix Diagonal. More...
 
static RCP< Vector > GetMatrixDiagonalInverse (const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool doLumped=false)
 Extract Matrix Diagonal. More...
 
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal (Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::magnitude(Teuchos::ScalarTraits< Scalar >::zero()), Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
 Extract Matrix Diagonal of lumped matrix. More...
 
static Teuchos::ArrayRCP
< Magnitude
GetMatrixMaxMinusOffDiagonal (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
 Return vector containing: max_{i=k}(-a_ik), for each for i in the matrix. More...
 
static Teuchos::ArrayRCP
< Magnitude
GetMatrixMaxMinusOffDiagonal (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber)
 
static Teuchos::RCP< Vector > GetInverse (Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
 Return vector containing inverse of input vector. More...
 
static RCP< Vector > GetMatrixOverlappedDiagonal (const Matrix &A)
 Extract Overlapped Matrix Diagonal. More...
 
static RCP< Vector > GetMatrixOverlappedDeletedRowsum (const Matrix &A)
 Extract Overlapped Matrix Deleted Rowsum. More...
 
static RCP< Xpetra::Vector
< Magnitude, LocalOrdinal,
GlobalOrdinal, Node > > 
GetMatrixOverlappedAbsDeletedRowsum (const Matrix &A)
 
static Teuchos::Array< MagnitudeResidualNorm (const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
 
static Teuchos::Array< MagnitudeResidualNorm (const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
 
static RCP< MultiVector > Residual (const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
 
static void Residual (const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS, MultiVector &Resid)
 
static Scalar PowerMethod (const Matrix &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
 Power method. More...
 
static Scalar PowerMethod (const Matrix &A, const RCP< Vector > &diagInvVec, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
 Power method. More...
 
static RCP< Teuchos::FancyOStreamMakeFancy (std::ostream &os)
 
static Teuchos::ScalarTraits
< Scalar >::magnitudeType 
Distance2 (const Teuchos::Array< Teuchos::ArrayRCP< const Scalar >> &v, LocalOrdinal i0, LocalOrdinal i1)
 Squared distance between two rows in a multivector. More...
 
static Teuchos::ScalarTraits
< Scalar >::magnitudeType 
Distance2 (const Teuchos::ArrayView< double > &weight, const Teuchos::Array< Teuchos::ArrayRCP< const Scalar >> &v, LocalOrdinal i0, LocalOrdinal i1)
 Weighted squared distance between two rows in a multivector. More...
 
static Teuchos::ArrayRCP
< const bool > 
DetectDirichletRows (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
 Detect Dirichlet rows. More...
 
static Kokkos::View< bool
*, typename
NO::device_type::memory_space > 
DetectDirichletRows_kokkos (const Matrix &A, const Magnitude &tol=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< SC >::magnitudeType >::zero(), const bool count_twos_as_dirichlet=false)
 Detect Dirichlet rows. More...
 
static Kokkos::View< bool
*, typename Kokkos::HostSpace > 
DetectDirichletRows_kokkos_host (const Matrix &A, const Magnitude &tol=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< SC >::magnitudeType >::zero(), const bool count_twos_as_dirichlet=false)
 
static Teuchos::ArrayRCP
< const bool > 
DetectDirichletRowsExt (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool &bHasZeroDiagonal, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
 Detect Dirichlet rows (extended version) More...
 
static void FindNonZeros (const Teuchos::ArrayRCP< const Scalar > vals, Teuchos::ArrayRCP< bool > nonzeros)
 Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon. More...
 
static void FindNonZeros (const typename Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type::t_dev_const_um vals, Kokkos::View< bool *, typename Node::device_type > nonzeros)
 Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon. More...
 
static void DetectDirichletColsAndDomains (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
 Detects Dirichlet columns & domains from a list of Dirichlet rows. More...
 
static void DetectDirichletColsAndDomains (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Kokkos::View< bool *, typename Node::device_type > &dirichletRows, Kokkos::View< bool *, typename Node::device_type > dirichletCols, Kokkos::View< bool *, typename Node::device_type > dirichletDomain)
 Detects Dirichlet columns & domains from a list of Dirichlet rows. More...
 
static void ApplyRowSumCriterion (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
 Apply Rowsum Criterion. More...
 
static void ApplyRowSumCriterion (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
 
static void ApplyRowSumCriterion (const Matrix &A, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, typename NO::device_type::memory_space > &dirichletRows)
 
static void ApplyRowSumCriterionHost (const Matrix &A, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, Kokkos::HostSpace > &dirichletRows)
 
static void ApplyRowSumCriterion (const Matrix &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, typename NO::device_type::memory_space > &dirichletRows)
 
static void ApplyRowSumCriterionHost (const Matrix &A, const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &BlockNumber, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType rowSumTol, Kokkos::View< bool *, Kokkos::HostSpace > &dirichletRows)
 
static Teuchos::ArrayRCP
< const bool > 
DetectDirichletCols (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
 Detect Dirichlet columns based on Dirichlet rows. More...
 
static Kokkos::View< bool
*, typename NO::device_type > 
DetectDirichletCols (const Matrix &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows)
 Detect Dirichlet columns based on Dirichlet rows. More...
 
static Scalar Frobenius (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B)
 Frobenius inner product of two matrices. More...
 
static void SetRandomSeed (const Teuchos::Comm< int > &comm)
 Set seed for random number generator. More...
 
static void FindDirichletRows (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, std::vector< LocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
 
static void ApplyOAZToMatrixRows (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, const std::vector< LocalOrdinal > &dirichletRows)
 
static void ApplyOAZToMatrixRows (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, const Teuchos::ArrayRCP< const bool > &dirichletRows)
 
static void ApplyOAZToMatrixRows (RCP< Matrix > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows)
 
static void ZeroDirichletRows (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, const std::vector< LocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
 
static void ZeroDirichletRows (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
 
static void ZeroDirichletRows (Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &X, const Teuchos::ArrayRCP< const bool > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
 
static void ZeroDirichletRows (RCP< Matrix > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
 
static void ZeroDirichletRows (RCP< MultiVector > &X, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
 
static void ZeroDirichletCols (Teuchos::RCP< Matrix > &A, const Teuchos::ArrayRCP< const bool > &dirichletCols, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
 
static void ZeroDirichletCols (RCP< Matrix > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletCols, SC replaceWith=Teuchos::ScalarTraits< SC >::zero())
 
static void FindDirichletRowsAndPropagateToCols (Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &A, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node >> &isDirichletRow, Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node >> &isDirichletCol)
 
static RCP< Xpetra::Matrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
ReplaceNonZerosWithOnes (const RCP< Matrix > &original)
 Creates a copy of a matrix where the non-zero entries are replaced by ones. More...
 
static RCP< const
Xpetra::BlockedMap
< LocalOrdinal, GlobalOrdinal,
Node > > 
GeneratedBlockedTargetMap (const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &sourceBlockedMap, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer)
 
static bool MapsAreNested (const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &colMap)
 
static RCP< Xpetra::Vector
< LocalOrdinal, LocalOrdinal,
GlobalOrdinal, Node > > 
CuthillMcKee (const Matrix &Op)
 
static RCP< Xpetra::Vector
< LocalOrdinal, LocalOrdinal,
GlobalOrdinal, Node > > 
ReverseCuthillMcKee (const Matrix &Op)
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >

MueLu utility class.

This class provides a number of static helper methods. Some are temporary and will eventually go away, while others should be moved to Xpetra.

Definition at line 51 of file MueLu_UtilitiesBase_fwd.hpp.

Member Typedef Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Magnitude

Definition at line 103 of file MueLu_UtilitiesBase_decl.hpp.

Member Function Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Crs2Op ( RCP< CrsMatrix >  Op)
static

Definition at line 80 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::CrsMatrixWrap< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetThresholdedMatrix ( const RCP< Matrix > &  Ain,
const Scalar  threshold,
const bool  keepDiagonal = true,
const GlobalOrdinal  expectedNNZperRow = -1 
)
static

Threshold a matrix.

Returns matrix filtered with a threshold value.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 89 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetThresholdedGraph ( const RCP< Matrix > &  A,
const Magnitude  threshold,
const GlobalOrdinal  expectedNNZperRow = -1 
)
static

Threshold a graph.

Returns graph filtered with a threshold value.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 138 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< Scalar > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixDiagonal_arcp ( const Matrix &  A)
static

Extract Matrix Diagonal.

Returns Matrix diagonal in ArrayRCP.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 171 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixDiagonal ( const Matrix &  A)
static

Extract Matrix Diagonal.

Returns Matrix diagonal in RCP<Vector>.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 196 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixDiagonalInverse ( const Matrix &  A,
Magnitude  tol = Teuchos::ScalarTraits<Scalar>::eps() * 100,
Scalar  valReplacement = Teuchos::ScalarTraits<Scalar>::zero(),
const bool  doLumped = false 
)
static

Extract Matrix Diagonal.

Returns inverse of the Matrix diagonal in ArrayRCP.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 233 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetLumpedMatrixDiagonal ( Matrix const &  A,
const bool  doReciprocal = false,
Magnitude  tol = Teuchos::ScalarTraits<Scalar>::magnitude(Teuchos::ScalarTraits<Scalar>::zero()),
Scalar  valReplacement = Teuchos::ScalarTraits<Scalar>::zero(),
const bool  replaceSingleEntryRowWithZero = false,
const bool  useAverageAbsDiagVal = false 
)
static

Extract Matrix Diagonal of lumped matrix.

Returns Matrix diagonal of lumped matrix in RCP<Vector>.

NOTE – it's assumed that A has been fillComplete'd.

Definition at line 324 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixMaxMinusOffDiagonal ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A)
static

Return vector containing: max_{i=k}(-a_ik), for each for i in the matrix.

Parameters
[in]A,:input matrix : vector containing max_{i=k}(-a_ik)

Definition at line 511 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixMaxMinusOffDiagonal ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  BlockNumber 
)
static

Definition at line 533 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetInverse ( Teuchos::RCP< const Vector >  v,
Magnitude  tol = Teuchos::ScalarTraits<Scalar>::eps() * 100,
Scalar  valReplacement = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Return vector containing inverse of input vector.

Parameters
[in]v,:input vector
[in]tol,:tolerance. If entries of input vector are smaller than tolerance they are replaced by valReplacement (see below). The default value for tol is 100*eps (machine precision)
[in]valReplacement,:Value put in for undefined entries in output vector (default: 0.0) : vector containing inverse values of input vector v

Definition at line 561 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixOverlappedDiagonal ( const Matrix &  A)
static

Extract Overlapped Matrix Diagonal.

Returns overlapped Matrix diagonal in ArrayRCP.

The local overlapped diagonal has an entry for each index in A's column map. NOTE – it's assumed that A has been fillComplete'd.

Definition at line 632 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixOverlappedDeletedRowsum ( const Matrix &  A)
static

Extract Overlapped Matrix Deleted Rowsum.

Returns overlapped Matrix deleted Rowsum in ArrayRCP.

The local overlapped deleted Rowsum has an entry for each index in A's column map. NOTE – it's assumed that A has been fillComplete'd.

Definition at line 648 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::Vector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetMatrixOverlappedAbsDeletedRowsum ( const Matrix &  A)
static

Definition at line 688 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::Array< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ResidualNorm ( const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Op,
const MultiVector &  X,
const MultiVector &  RHS 
)
static

Definition at line 731 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::Array< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ResidualNorm ( const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Op,
const MultiVector &  X,
const MultiVector &  RHS,
MultiVector &  Resid 
)
static

Definition at line 743 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Residual ( const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Op,
const MultiVector &  X,
const MultiVector &  RHS 
)
static

Definition at line 755 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Residual ( const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Op,
const MultiVector &  X,
const MultiVector &  RHS,
MultiVector &  Resid 
)
static

Definition at line 766 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static Scalar MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PowerMethod ( const Matrix &  A,
bool  scaleByDiag = true,
LocalOrdinal  niters = 10,
Magnitude  tolerance = 1e-2,
bool  verbose = false,
unsigned int  seed = 123 
)
static

Power method.

Parameters
Amatrix
scaleByDiagif true, estimate the largest eigenvalue of \( D^; A \).
nitersmaximum number of iterations
tolerancestopping tolerance if true, print iteration information seed for randomizing initial guess

(Shamelessly grabbed from tpetra/examples.)

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static Scalar MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::PowerMethod ( const Matrix &  A,
const RCP< Vector > &  diagInvVec,
LocalOrdinal  niters = 10,
Magnitude  tolerance = 1e-2,
bool  verbose = false,
unsigned int  seed = 123 
)
static

Power method.

Parameters
Amatrix
diagInvVecreciprocal of matrix diagonal
nitersmaximum number of iterations
tolerancestopping tolerance if true, print iteration information seed for randomizing initial guess

(Shamelessly grabbed from tpetra/examples.)

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Teuchos::FancyOStream > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::MakeFancy ( std::ostream &  os)
static

Definition at line 844 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
Teuchos::ScalarTraits< Scalar >::magnitudeType MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Distance2 ( const Teuchos::Array< Teuchos::ArrayRCP< const Scalar >> &  v,
LocalOrdinal  i0,
LocalOrdinal  i1 
)
static

Squared distance between two rows in a multivector.

Used for coordinate vectors.

Definition at line 852 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
Teuchos::ScalarTraits< Scalar >::magnitudeType MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Distance2 ( const Teuchos::ArrayView< double > &  weight,
const Teuchos::Array< Teuchos::ArrayRCP< const Scalar >> &  v,
LocalOrdinal  i0,
LocalOrdinal  i1 
)
static

Weighted squared distance between two rows in a multivector.

Used for coordinate vectors.

Definition at line 865 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< const bool > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletRows ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Magnitude tol = Teuchos::ScalarTraits<Magnitude>::zero(),
bool  count_twos_as_dirichlet = false 
)
static

Detect Dirichlet rows.

The routine assumes, that if there is only one nonzero per row, it is on the diagonal and therefore a DBC. This is safe for most of our applications, but one should be aware of that.

There is an alternative routine (see DetectDirichletRowsExt)

Parameters
[in]Amatrix
[in]tolIf a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero.
Returns
boolean array. The ith entry is true iff row i is a Dirichlet row.

Definition at line 879 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Kokkos::View< bool *, typename NO::device_type::memory_space > MueLu::UtilitiesBase< SC, LO, GO, NO >::DetectDirichletRows_kokkos ( const Matrix &  A,
const Magnitude tol = Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<SC>::magnitudeType>::zero(),
const bool  count_twos_as_dirichlet = false 
)
static

Detect Dirichlet rows.

Parameters
[in]Amatrix
[in]tolIf a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero.
Returns
boolean array. The ith entry is true iff row i is a Dirichlet row.

Definition at line 1023 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Kokkos::View< bool *, typename Kokkos::HostSpace > MueLu::UtilitiesBase< SC, LO, GO, NO >::DetectDirichletRows_kokkos_host ( const Matrix &  A,
const Magnitude tol = Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<SC>::magnitudeType>::zero(),
const bool  count_twos_as_dirichlet = false 
)
static

Definition at line 1032 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< const bool > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletRowsExt ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
bool &  bHasZeroDiagonal,
const Magnitude tol = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Detect Dirichlet rows (extended version)

Look at each matrix row and mark it as Dirichlet if there is only one "not small" nonzero on the diagonal. In determining whether a nonzero is "not small" use abs(A(i,j)) / sqrt(abs(diag[i]*diag[j])) > tol

Parameters
[in]Amatrix
in/out]bHasZeroDiagonal Reference to boolean variable. Returns true if there is a zero on the diagonal in the local part of the Matrix. Otherwise it is false. Different processors might return a different value. There is no global reduction!
[in]tolIf a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero.
Returns
boolean array. The ith entry is true iff row i is a Dirichlet row.

Definition at line 1041 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FindNonZeros ( const Teuchos::ArrayRCP< const Scalar vals,
Teuchos::ArrayRCP< bool >  nonzeros 
)
static

Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.

Parameters
[in]vals- ArrayRCP<const Scalar> of values to be tested
[out]nonzeros- ArrayRCP<bool> of true/false values for whether each entry in vals is nonzero

Definition at line 1076 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FindNonZeros ( const typename Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type::t_dev_const_um  vals,
Kokkos::View< bool *, typename Node::device_type >  nonzeros 
)
static

Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.

Parameters
[in]vals- ArrayRCP<const Scalar> of values to be tested
[out]nonzeros- ArrayRCP<bool> of true/false values for whether each entry in vals is nonzero

Definition at line 1089 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletColsAndDomains ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Teuchos::ArrayRCP< bool > &  dirichletRows,
Teuchos::ArrayRCP< bool >  dirichletCols,
Teuchos::ArrayRCP< bool >  dirichletDomain 
)
static

Detects Dirichlet columns & domains from a list of Dirichlet rows.

Parameters
[in]A- Matrix on which to apply Dirichlet column detection
[in]dirichletRows- ArrayRCP<bool> of indicators as to which rows are Dirichlet
[out]dirichletCols- ArrayRCP<bool> of indicators as to which cols are Dirichlet
[out]dirichletDomain- ArrayRCP<bool> of indicators as to which domains are Dirichlet

Definition at line 1106 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletColsAndDomains ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Kokkos::View< bool *, typename Node::device_type > &  dirichletRows,
Kokkos::View< bool *, typename Node::device_type >  dirichletCols,
Kokkos::View< bool *, typename Node::device_type >  dirichletDomain 
)
static

Detects Dirichlet columns & domains from a list of Dirichlet rows.

Parameters
[in]A- Matrix on which to apply Dirichlet column detection
[in]dirichletRows- View<bool> of indicators as to which rows are Dirichlet
[out]dirichletCols- View<bool> of indicators as to which cols are Dirichlet
[out]dirichletDomain- View<bool> of indicators as to which domains are Dirichlet

Definition at line 1146 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterion ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Magnitude  rowSumTol,
Teuchos::ArrayRCP< bool > &  dirichletRows 
)
static

Apply Rowsum Criterion.

Flags a row i as dirichlet if:

{j=i} A_ij > A_ii * tol

Parameters
[in]Amatrix
[in]rowSumTolSee above
in/out]dirichletRows boolean array. The ith entry is true if the above criterion is satisfied (or if it was already set to true)
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterion ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  BlockNumber,
const Magnitude  rowSumTol,
Teuchos::ArrayRCP< bool > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterion ( const Matrix &  A,
const typename Teuchos::ScalarTraits< Scalar >::magnitudeType  rowSumTol,
Kokkos::View< bool *, typename NO::device_type::memory_space > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterionHost ( const Matrix &  A,
const typename Teuchos::ScalarTraits< Scalar >::magnitudeType  rowSumTol,
Kokkos::View< bool *, Kokkos::HostSpace > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterion ( const Matrix &  A,
const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  BlockNumber,
const typename Teuchos::ScalarTraits< Scalar >::magnitudeType  rowSumTol,
Kokkos::View< bool *, typename NO::device_type::memory_space > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyRowSumCriterionHost ( const Matrix &  A,
const Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  BlockNumber,
const typename Teuchos::ScalarTraits< Scalar >::magnitudeType  rowSumTol,
Kokkos::View< bool *, Kokkos::HostSpace > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ArrayRCP< const bool > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletCols ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Teuchos::ArrayRCP< const bool > &  dirichletRows 
)
static

Detect Dirichlet columns based on Dirichlet rows.

The routine finds all column indices that are in Dirichlet rows, where Dirichlet rows are described by dirichletRows, as returned by DetectDirichletRows.

Parameters
[in]Amatrix
[in]dirichletRowsarray of Dirichlet rows.
Returns
boolean array. The ith entry is true iff row i is a Dirichlet column.

Definition at line 1360 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static Kokkos::View<bool*, typename NO::device_type> MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DetectDirichletCols ( const Matrix &  A,
const Kokkos::View< const bool *, typename NO::device_type > &  dirichletRows 
)
static

Detect Dirichlet columns based on Dirichlet rows.

The routine finds all column indices that are in Dirichlet rows, where Dirichlet rows are described by dirichletRows, as returned by DetectDirichletRows.

Parameters
[in]Amatrix
[in]dirichletRowsarray of Dirichlet rows.
Returns
boolean array. The ith entry is true iff row i is a Dirichlet column.
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Scalar MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Frobenius ( const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B 
)
static

Frobenius inner product of two matrices.

Used in energy minimization algorithms

Definition at line 1451 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetRandomSeed ( const Teuchos::Comm< int > &  comm)
static

Set seed for random number generator.

Distribute the seeds evenly in [1,INT_MAX-1]. This guarantees nothing about where random number streams on difference processes will intersect. This does avoid overflow situations in parallel when multiplying by a PID. It also avoids the pathological case of having the same random number stream on each process.

Definition at line 1511 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FindDirichletRows ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
std::vector< LocalOrdinal > &  dirichletRows,
bool  count_twos_as_dirichlet = false 
)
static

Definition at line 1535 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyOAZToMatrixRows ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
const std::vector< LocalOrdinal > &  dirichletRows 
)
static

Definition at line 1557 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyOAZToMatrixRows ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
const Teuchos::ArrayRCP< const bool > &  dirichletRows 
)
static

Definition at line 1583 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ApplyOAZToMatrixRows ( RCP< Matrix > &  A,
const Kokkos::View< const bool *, typename Node::device_type > &  dirichletRows 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletRows ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
const std::vector< LocalOrdinal > &  dirichletRows,
Scalar  replaceWith = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Definition at line 1655 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletRows ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
const Teuchos::ArrayRCP< const bool > &  dirichletRows,
Scalar  replaceWith = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Definition at line 1671 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletRows ( Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  X,
const Teuchos::ArrayRCP< const bool > &  dirichletRows,
Scalar  replaceWith = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Definition at line 1690 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletRows ( RCP< Matrix > &  A,
const Kokkos::View< const bool *, typename NO::device_type > &  dirichletRows,
SC  replaceWith = Teuchos::ScalarTraitsSC >::zero() 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletRows ( RCP< MultiVector > &  X,
const Kokkos::View< const bool *, typename NO::device_type > &  dirichletRows,
SC  replaceWith = Teuchos::ScalarTraitsSC >::zero() 
)
static
template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletCols ( Teuchos::RCP< Matrix > &  A,
const Teuchos::ArrayRCP< const bool > &  dirichletCols,
Scalar  replaceWith = Teuchos::ScalarTraits<Scalar>::zero() 
)
static

Definition at line 1751 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
static void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ZeroDirichletCols ( RCP< Matrix > &  A,
const Kokkos::View< const bool *, typename NO::device_type > &  dirichletCols,
SC  replaceWith = Teuchos::ScalarTraitsSC >::zero() 
)
static
template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FindDirichletRowsAndPropagateToCols ( Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &  A,
Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node >> &  isDirichletRow,
Teuchos::RCP< Xpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node >> &  isDirichletCol 
)
static

Definition at line 1794 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReplaceNonZerosWithOnes ( const RCP< Matrix > &  original)
static

Creates a copy of a matrix where the non-zero entries are replaced by ones.

Definition at line 1837 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GeneratedBlockedTargetMap ( const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &  sourceBlockedMap,
const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &  Importer 
)
static

Definition at line 1869 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
bool MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::MapsAreNested ( const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &  rowMap,
const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &  colMap 
)
static

Definition at line 1921 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CuthillMcKee ( const Matrix &  Op)
static

Perform a Cuthill-McKee (CM) or Reverse Cuthill-McKee (RCM) ordering of the local component of the matrix. Kokkos-Kernels has an RCM implementation, so we reverse that here if we call CM.

Definition at line 1972 of file MueLu_UtilitiesBase_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReverseCuthillMcKee ( const Matrix &  Op)
static

Perform a Reverse Cuthill-McKee (RCM) ordering of the local component of the matrix.

Definition at line 1943 of file MueLu_UtilitiesBase_def.hpp.


The documentation for this class was generated from the following files: