Belos Package Browser (Single Doxygen Collection)
Development
|
Simple example of a user's defined Belos::MultiVec class. More...
#include <MyMultiVec.hpp>
Public Member Functions | |
MyMultiVec (const ptrdiff_t Length, const int NumberVecs) | |
Constructor for a NumberVecs vectors of length Length . More... | |
MyMultiVec (const ptrdiff_t Length, const std::vector< ScalarType * > &rhs) | |
Constructor with already allocated memory. More... | |
MyMultiVec (const MyMultiVec &rhs) | |
Copy constructor, performs a deep copy. More... | |
~MyMultiVec () | |
Destructor. More... | |
MyMultiVec * | Clone (const int NumberVecs) const |
Returns a clone of the current std::vector. More... | |
MyMultiVec * | CloneCopy () const |
Create a new MultiVec and copy contents of *this into it (deep copy). More... | |
MyMultiVec * | CloneCopy (const std::vector< int > &index) const |
Returns a clone copy of specified vectors. More... | |
MyMultiVec * | CloneViewNonConst (const std::vector< int > &index) |
Returns a view of current std::vector (shallow copy) More... | |
const MyMultiVec * | CloneView (const std::vector< int > &index) const |
Returns a view of current std::vector (shallow copy), const version. More... | |
ptrdiff_t | GetGlobalLength () const |
The number of rows in the multivector. More... | |
int | GetNumberVecs () const |
The number of vectors (i.e., columns) in the multivector. More... | |
void | MvTimesMatAddMv (const ScalarType alpha, const Belos::MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta) |
Update *this with alpha * A * B + beta * (*this ). More... | |
void | MvAddMv (const ScalarType alpha, const Belos::MultiVec< ScalarType > &A, const ScalarType beta, const Belos::MultiVec< ScalarType > &B) |
Replace *this with alpha * A + beta * B . More... | |
void | MvScale (const ScalarType alpha) |
Scale each element of the vectors in *this with alpha . More... | |
void | MvScale (const std::vector< ScalarType > &alpha) |
Scale each element of the i -th vector in *this with alpha[i] . More... | |
void | MvTransMv (const ScalarType alpha, const Belos::MultiVec< ScalarType > &A, Teuchos::SerialDenseMatrix< int, ScalarType > &B) const |
Compute a dense matrix B through the matrix-matrix multiply alpha * A^T * (*this ). More... | |
void | MvDot (const Belos::MultiVec< ScalarType > &A, std::vector< ScalarType > &b) const |
Compute the dot product of each column of *this with the corresponding column of A. More... | |
void | MvNorm (std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, Belos::NormType type=Belos::TwoNorm) const |
Compute the norm of each vector in *this . More... | |
void | SetBlock (const Belos::MultiVec< ScalarType > &A, const std::vector< int > &index) |
Copy the vectors in A to a set of vectors in *this . More... | |
void | MvRandom () |
Fill all the vectors in *this with random numbers. More... | |
void | MvInit (const ScalarType alpha) |
Replace each element of the vectors in *this with alpha . More... | |
void | MvPrint (std::ostream &os) const |
Print *this multivector to the os output stream. More... | |
ScalarType & | operator() (const int i, const int j) |
const ScalarType & | operator() (const int i, const int j) const |
ScalarType * | operator[] (int v) |
ScalarType * | operator[] (int v) const |
Public Member Functions inherited from Belos::MultiVec< ScalarType > | |
MultiVec () | |
Default constructor. More... | |
virtual | ~MultiVec () |
Destructor (virtual for memory safety of derived classes). More... | |
Private Member Functions | |
void | Check () |
Private Attributes | |
const ptrdiff_t | Length_ |
Length of the vectors. More... | |
const int | NumberVecs_ |
Number of multi-vectors. More... | |
std::vector< ScalarType * > | data_ |
Pointers to the storage of the vectors. More... | |
std::vector< bool > | ownership_ |
If true , then this object owns the vectors and must free them in dtor. More... | |
Simple example of a user's defined Belos::MultiVec class.
This is a simple, single processor example of user's defined MultiVec-derived class. The class is templated with ScalarType; possible choices are, for example, "float", "double", or "std::complex<double>".
Definition at line 65 of file MyMultiVec.hpp.
|
inline |
Constructor for a NumberVecs
vectors of length Length
.
Definition at line 70 of file MyMultiVec.hpp.
|
inline |
Constructor with already allocated memory.
Definition at line 90 of file MyMultiVec.hpp.
|
inline |
Copy constructor, performs a deep copy.
Definition at line 108 of file MyMultiVec.hpp.
|
inline |
Destructor.
Definition at line 131 of file MyMultiVec.hpp.
|
inlinevirtual |
Returns a clone of the current std::vector.
Implements Belos::MultiVec< ScalarType >.
Definition at line 139 of file MyMultiVec.hpp.
|
inlinevirtual |
Create a new MultiVec and copy contents of *this
into it (deep copy).
Implements Belos::MultiVec< ScalarType >.
Definition at line 152 of file MyMultiVec.hpp.
|
inlinevirtual |
Returns a clone copy of specified vectors.
Implements Belos::MultiVec< ScalarType >.
Definition at line 158 of file MyMultiVec.hpp.
|
inlinevirtual |
Returns a view of current std::vector (shallow copy)
Implements Belos::MultiVec< ScalarType >.
Definition at line 173 of file MyMultiVec.hpp.
|
inlinevirtual |
Returns a view of current std::vector (shallow copy), const version.
Implements Belos::MultiVec< ScalarType >.
Definition at line 185 of file MyMultiVec.hpp.
|
inlinevirtual |
The number of rows in the multivector.
Implements Belos::MultiVec< ScalarType >.
Definition at line 197 of file MyMultiVec.hpp.
|
inlinevirtual |
The number of vectors (i.e., columns) in the multivector.
Implements Belos::MultiVec< ScalarType >.
Definition at line 202 of file MyMultiVec.hpp.
|
inlinevirtual |
Update *this
with alpha
* A
* B
+ beta
* (*this
).
Implements Belos::MultiVec< ScalarType >.
Definition at line 208 of file MyMultiVec.hpp.
|
inlinevirtual |
Replace *this
with alpha
* A
+ beta
* B
.
Implements Belos::MultiVec< ScalarType >.
Definition at line 259 of file MyMultiVec.hpp.
|
inlinevirtual |
Scale each element of the vectors in *this
with alpha
.
Implements Belos::MultiVec< ScalarType >.
Definition at line 282 of file MyMultiVec.hpp.
|
inlinevirtual |
Scale each element of the i
-th vector in *this
with alpha[i]
.
Implements Belos::MultiVec< ScalarType >.
Definition at line 292 of file MyMultiVec.hpp.
|
inlinevirtual |
Compute a dense matrix B
through the matrix-matrix multiply alpha
* A^T
* (*this
).
Implements Belos::MultiVec< ScalarType >.
Definition at line 303 of file MyMultiVec.hpp.
|
inlinevirtual |
Compute the dot product of each column of *this with the corresponding column of A.
Compute a vector b
whose entries are the individual dot-products. That is, b[i] = A[i]^H * (*this)[i]
where A[i]
is the i-th column of A.
Implements Belos::MultiVec< ScalarType >.
Definition at line 327 of file MyMultiVec.hpp.
|
inlinevirtual |
Compute the norm of each vector in *this
.
normvec | [out] On output, normvec[i] holds the norm of the i-th vector of *this . |
type | [in] The type of norm to compute. The 2-norm is the default. |
Implements Belos::MultiVec< ScalarType >.
Definition at line 347 of file MyMultiVec.hpp.
|
inlinevirtual |
Copy the vectors in A
to a set of vectors in *this
.
The numvecs
vectors in A
are copied to a subset of vectors in *this
indicated by the indices given in index
.
Implements Belos::MultiVec< ScalarType >.
Definition at line 368 of file MyMultiVec.hpp.
|
inlinevirtual |
Fill all the vectors in *this
with random numbers.
Implements Belos::MultiVec< ScalarType >.
Definition at line 386 of file MyMultiVec.hpp.
|
inlinevirtual |
Replace each element of the vectors in *this
with alpha
.
Implements Belos::MultiVec< ScalarType >.
Definition at line 398 of file MyMultiVec.hpp.
|
inlinevirtual |
Print *this
multivector to the os
output stream.
Implements Belos::MultiVec< ScalarType >.
Definition at line 407 of file MyMultiVec.hpp.
|
inline |
Definition at line 421 of file MyMultiVec.hpp.
|
inline |
Definition at line 429 of file MyMultiVec.hpp.
|
inline |
Definition at line 437 of file MyMultiVec.hpp.
|
inline |
Definition at line 443 of file MyMultiVec.hpp.
|
inlineprivate |
Definition at line 449 of file MyMultiVec.hpp.
|
private |
Length of the vectors.
Definition at line 459 of file MyMultiVec.hpp.
|
private |
Number of multi-vectors.
Definition at line 461 of file MyMultiVec.hpp.
|
private |
Pointers to the storage of the vectors.
Definition at line 463 of file MyMultiVec.hpp.
|
private |
If true
, then this object owns the vectors and must free them in dtor.
Definition at line 465 of file MyMultiVec.hpp.