MueLu
Version of the Day
|
#include <MueLu_UtilitiesBase_fwd.hpp>
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 Magnitude 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< Scalar > | GetMatrixDiagonal_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::RCP< Vector > | 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::RCP< Vector > | 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 GlobalOrdinal | CountNegativeDiagonalEntries (const Matrix &A) |
Counts the number of negative diagonal entries. More... | |
static Teuchos::Array< Magnitude > | ResidualNorm (const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS) |
static Teuchos::Array< Magnitude > | ResidualNorm (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::FancyOStream > | MakeFancy (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 | EnforceInitialCondition (const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &InitialGuess, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), const bool count_twos_as_dirichlet=false) |
Detect Dirichlet rows and copy values from RHS multivector to InitialGuess for Dirichlet rows. 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) |
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 15 of file MueLu_UtilitiesBase_fwd.hpp.
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType MueLu::UtilitiesBase< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Magnitude |
Definition at line 71 of file MueLu_UtilitiesBase_decl.hpp.
|
static |
Definition at line 45 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Threshold a matrix.
Returns matrix filtered with a threshold value.
NOTE – it's assumed that A has been fillComplete'd.
Definition at line 169 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Threshold a graph.
Returns graph filtered with a threshold value.
NOTE – it's assumed that A has been fillComplete'd.
Definition at line 225 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Extract Matrix Diagonal.
Returns Matrix diagonal in ArrayRCP.
NOTE – it's assumed that A has been fillComplete'd.
Definition at line 258 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Extract Matrix Diagonal.
Returns Matrix diagonal in RCP<Vector>.
NOTE – it's assumed that A has been fillComplete'd.
Definition at line 283 of file MueLu_UtilitiesBase_def.hpp.
|
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 319 of file MueLu_UtilitiesBase_def.hpp.
|
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 414 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Return vector containing: max_{i=k}(-a_ik), for each for i in the matrix.
[in] | A,: | input matrix : vector containing max_{i=k}(-a_ik) |
Definition at line 613 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 654 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Return vector containing inverse of input vector.
[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 698 of file MueLu_UtilitiesBase_def.hpp.
|
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 769 of file MueLu_UtilitiesBase_def.hpp.
|
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 787 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 827 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Counts the number of negative diagonal entries.
Returns a GlobalOrdinal with the number of negative diagonal entries This generally will involve MPI communication and this must be called on all ranks in A's communicator. NOTE: This only works on matrices locally fitted column maps.
Definition at line 870 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 901 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 913 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 925 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 936 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Power method.
A | matrix |
scaleByDiag | if true, estimate the largest eigenvalue of \( D^; A \). |
niters | maximum number of iterations |
tolerance | stopping tolerance if true, print iteration information seed for randomizing initial guess |
(Shamelessly grabbed from tpetra/examples.)
|
static |
Power method.
A | matrix |
diagInvVec | reciprocal of matrix diagonal |
niters | maximum number of iterations |
tolerance | stopping tolerance if true, print iteration information seed for randomizing initial guess |
(Shamelessly grabbed from tpetra/examples.)
|
static |
Definition at line 1014 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Squared distance between two rows in a multivector.
Used for coordinate vectors.
Definition at line 1022 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Weighted squared distance between two rows in a multivector.
Used for coordinate vectors.
Definition at line 1035 of file MueLu_UtilitiesBase_def.hpp.
|
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)
[in] | A | matrix |
[in] | tol | If a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero. |
Definition at line 1049 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Detect Dirichlet rows.
[in] | A | matrix |
[in] | tol | If a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero. |
Definition at line 1213 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1222 of file MueLu_UtilitiesBase_def.hpp.
|
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
[in] | A | matrix |
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] | tol | If a row entry's magnitude is less than or equal to this tolerance, the entry is treated as zero. |
Definition at line 1231 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Detect Dirichlet rows and copy values from RHS multivector to InitialGuess for Dirichlet rows.
This can be used to assure that the InitialGuess satisfies the boundary conditions enforced on A. Useful in particular for using CG when boundary conditions have only been enforce by one-and-zeroing rows of A, but not columns.
Definition at line 1266 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.
[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 1297 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Find non-zero values in an ArrayRCP Compares the value to 2 * machine epsilon.
[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 1310 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Detects Dirichlet columns & domains from a list of Dirichlet rows.
[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 1335 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Detects Dirichlet columns & domains from a list of Dirichlet rows.
[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 1375 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Apply Rowsum Criterion.
Flags a row i as dirichlet if:
{j=i} A_ij > A_ii * tol
[in] | A | matrix |
[in] | rowSumTol | See 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) |
|
static |
|
static |
|
static |
|
static |
|
static |
|
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.
[in] | A | matrix |
[in] | dirichletRows | array of Dirichlet rows. |
Definition at line 1597 of file MueLu_UtilitiesBase_def.hpp.
|
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.
[in] | A | matrix |
[in] | dirichletRows | array of Dirichlet rows. |
|
static |
Frobenius inner product of two matrices.
Used in energy minimization algorithms
Definition at line 1696 of file MueLu_UtilitiesBase_def.hpp.
|
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 1756 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1780 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1802 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1828 of file MueLu_UtilitiesBase_def.hpp.
|
static |
|
static |
Definition at line 1908 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1924 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 1943 of file MueLu_UtilitiesBase_def.hpp.
|
static |
|
static |
|
static |
Definition at line 2012 of file MueLu_UtilitiesBase_def.hpp.
|
static |
|
static |
Definition at line 2059 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Creates a copy of a matrix where the non-zero entries are replaced by ones.
Definition at line 2102 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 2146 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Definition at line 2198 of file MueLu_UtilitiesBase_def.hpp.
|
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 2255 of file MueLu_UtilitiesBase_def.hpp.
|
static |
Perform a Reverse Cuthill-McKee (RCM) ordering of the local component of the matrix.
Definition at line 2226 of file MueLu_UtilitiesBase_def.hpp.