Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
|
Partial specialization of MultiVecTraits for MV = Tpetra::MultiVector. More...
#include <Belos_TpetraAdapter_MP_Vector.hpp>
Public Types | |
typedef Storage::ordinal_type | s_ordinal |
typedef Storage::value_type | BaseScalar |
typedef Sacado::MP::Vector< Storage > | Scalar |
typedef Tpetra::MultiVector< Scalar, LO, GO, Node > | MV |
typedef Tpetra::MultiVector< Scalar, LO, GO, Node >::dot_type | dot_type |
typedef Tpetra::MultiVector< Scalar, LO, GO, Node >::mag_type | mag_type |
Static Public Member Functions | |
static Teuchos::RCP< MV > | Clone (const MV &X, const int numVecs) |
Create a new MultiVector with numVecs columns. More... | |
static Teuchos::RCP< MV > | CloneCopy (const MV &X) |
Create and return a deep copy of X. More... | |
static Teuchos::RCP< MV > | CloneCopy (const MV &mv, const std::vector< int > &index) |
Create and return a deep copy of the given columns of mv. More... | |
static Teuchos::RCP< MV > | CloneCopy (const MV &mv, const Teuchos::Range1D &index) |
Create and return a deep copy of the given columns of mv. More... | |
static Teuchos::RCP< MV > | CloneViewNonConst (MV &mv, const std::vector< int > &index) |
static Teuchos::RCP< MV > | CloneViewNonConst (MV &mv, const Teuchos::Range1D &index) |
static Teuchos::RCP< const MV > | CloneView (const MV &mv, const std::vector< int > &index) |
static Teuchos::RCP< const MV > | CloneView (const MV &mv, const Teuchos::Range1D &index) |
static ptrdiff_t | GetGlobalLength (const MV &mv) |
static int | GetNumberVecs (const MV &mv) |
static bool | HasConstantStride (const MV &mv) |
static void | MvTimesMatAddMv (const dot_type &alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Teuchos::SerialDenseMatrix< int, dot_type > &B, const dot_type &beta, Tpetra::MultiVector< Scalar, LO, GO, Node > &C) |
static void | MvAddMv (Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, Scalar beta, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv) |
mv := alpha*A + beta*B More... | |
static void | MvScale (Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Scalar &alpha) |
static void | MvScale (Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< BaseScalar > &alphas) |
static void | MvScale (Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< Scalar > &alphas) |
static void | MvTransMv (dot_type alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Teuchos::SerialDenseMatrix< int, dot_type > &C) |
static void | MvDot (const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, std::vector< dot_type > &dots) |
For all columns j of A, set dots[j] := A[j]^T * B[j] . More... | |
static void | MvNorm (const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, std::vector< mag_type > &normvec, NormType type=TwoNorm) |
For all columns j of mv, set normvec[j] = norm(mv[j]) . More... | |
static void | SetBlock (const MV &A, const std::vector< int > &index, MV &mv) |
static void | SetBlock (const MV &A, const Teuchos::Range1D &index, MV &mv) |
static void | Assign (const MV &A, MV &mv) |
static void | MvRandom (MV &mv) |
static void | MvInit (MV &mv, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::zero()) |
static void | MvPrint (const MV &mv, std::ostream &os) |
Partial specialization of MultiVecTraits for MV = Tpetra::MultiVector.
This interface lets Belos' solvers work directly with Tpetra::MultiVector objects as the multivector type (corresponding to the MV template parameter).
The four template parameters of this partial specialization correspond exactly to the four template parameters of Tpetra::MultiVector. See the Tpetra::MultiVector documentation for more information.
Definition at line 73 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Storage::ordinal_type Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::s_ordinal |
Definition at line 77 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Storage::value_type Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::BaseScalar |
Definition at line 78 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Sacado::MP::Vector<Storage> Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::Scalar |
Definition at line 79 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Tpetra::MultiVector<Scalar, LO, GO, Node> Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::MV |
Definition at line 80 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Tpetra::MultiVector<Scalar,LO,GO,Node>::dot_type Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::dot_type |
Definition at line 82 of file Belos_TpetraAdapter_MP_Vector.hpp.
typedef Tpetra::MultiVector<Scalar,LO,GO,Node>::mag_type Belos::MultiVecTraits< typename Storage::value_type, Tpetra::MultiVector< Sacado::MP::Vector< Storage >, LO, GO, Node > >::mag_type |
Definition at line 83 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Create a new MultiVector with numVecs
columns.
The returned Tpetra::MultiVector has the same Tpetra::Map (distribution over one or more parallel processes) as X
. Its entries are not initialized and have undefined values.
Definition at line 90 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Create and return a deep copy of X.
Definition at line 97 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Create and return a deep copy of the given columns of mv.
0 <= k < index.size()
, Definition at line 125 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Create and return a deep copy of the given columns of mv.
Definition at line 161 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 191 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 222 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 255 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 286 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 315 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 319 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 323 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 328 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
mv := alpha*A + beta*B
The Tpetra specialization of this method ignores and completely overwrites any NaN or Inf entries in A. Thus, it does not mean the same thing as mv := 0*mv + alpha*A + beta*B
in IEEE 754 floating-point arithmetic. (Remember that NaN*0 = NaN.)
Definition at line 396 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 406 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 413 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 424 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 431 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
For all columns j of A, set dots[j] := A[j]^T * B[j]
.
Definition at line 519 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
For all columns j of mv, set normvec[j] = norm(mv[j])
.
Definition at line 546 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 589 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 612 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 687 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 734 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 739 of file Belos_TpetraAdapter_MP_Vector.hpp.
|
inlinestatic |
Definition at line 744 of file Belos_TpetraAdapter_MP_Vector.hpp.