Tpetra parallel linear algebra
Version of the Day
|
A distributed dense vector. More...
#include <Tpetra_Vector_decl.hpp>
Public Types | |
Typedefs to facilitate template metaprogramming | |
typedef Scalar | scalar_type |
This class' first template parameter; the type of each entry in the Vector. More... | |
typedef base_type::impl_scalar_type | impl_scalar_type |
The type used internally in place of Scalar . More... | |
typedef LocalOrdinal | local_ordinal_type |
This class' second template parameter; the type of local indices. More... | |
typedef GlobalOrdinal | global_ordinal_type |
This class' third template parameter; the type of global indices. More... | |
typedef Node::device_type | device_type |
The Kokkos device type. More... | |
typedef Node | node_type |
The Kokkos Node type. More... | |
typedef base_type::dot_type | dot_type |
Type of an inner ("dot") product result. More... | |
typedef base_type::mag_type | mag_type |
Type of a norm result. More... | |
typedef base_type::dual_view_type | dual_view_type |
Kokkos::DualView specialization used by this class. More... | |
typedef base_type::wrapped_dual_view_type | wrapped_dual_view_type |
WrappedDualView specialization used by this class. More... | |
typedef base_type::map_type | map_type |
The type of the Map specialization used by this class. More... | |
Typedefs to facilitate template metaprogramming. | |
using | execution_space = typename device_type::execution_space |
Type of the (new) Kokkos execution space. More... | |
using | host_view_type = typename dual_view_type::t_host |
using | device_view_type = typename dual_view_type::t_dev |
Typedefs | |
using | packet_type = typename::Kokkos::ArithTraits< Scalar >::val_type |
The type of each datum being sent or received in an Import or Export. More... | |
Public Member Functions | |
virtual void | removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) override |
Remove processes owning zero rows from the Map and their communicator. More... | |
void | setCopyOrView (const Teuchos::DataAccess copyOrView) |
Set whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics. More... | |
Teuchos::DataAccess | getCopyOrView () const |
Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics. More... | |
void | assign (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &src) |
Copy the contents of src into *this (deep copy). More... | |
Teuchos::RCP< MultiVector< T, LocalOrdinal, GlobalOrdinal, Node > > | convert () const |
Return another MultiVector with the same entries, but converted to a different Scalar type T . More... | |
bool | isSameSize (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const |
Constructors and destructor | |
Vector () | |
Default constructor: makes a Vector with no rows or columns. More... | |
Vector (const Teuchos::RCP< const map_type > &map, const bool zeroOut=true) | |
Basic constructor. More... | |
Vector (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::DataAccess copyOrView) | |
Copy constructor (shallow or deep copy). More... | |
Vector (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::RCP< const map_type > &map, const local_ordinal_type rowOffset=0) | |
"Offset view" constructor that views the input Vector's local data, but with the given Map, using the given row offset. More... | |
Vector (const Teuchos::RCP< const map_type > &map, const Teuchos::ArrayView< const Scalar > &A) | |
Set vector values from an existing array (copy) More... | |
Vector (const Teuchos::RCP< const map_type > &map, const dual_view_type &view) | |
Expert mode constructor, that takes a Kokkos::DualView of the Vector's data, and returns a Vector that views those data. More... | |
Vector (const Teuchos::RCP< const map_type > &map, const dual_view_type &view, const dual_view_type &origView) | |
Expert mode constructor, that takes a Kokkos::DualView of the Vector's data and the "original" Kokkos::DualView of the data, and returns a Vector that views those data. More... | |
Vector (const Teuchos::RCP< const map_type > &map, const wrapped_dual_view_type &d_view) | |
Expert mode constructor, that takes a WrappedDualView of the MultiVector's data. More... | |
Vector (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const size_t j) | |
Create a Vector that views a single column of the input MultiVector. More... | |
Vector (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default | |
Copy constructor (shallow copy). More... | |
Vector (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default | |
Move constructor (shallow move). More... | |
Vector & | operator= (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default |
Copy assignment (shallow copy). More... | |
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | operator= (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default |
Move assigment (shallow move). More... | |
virtual | ~Vector ()=default |
Destructor (virtual for memory safety of derived classes). More... | |
Post-construction modification routines | |
void | replaceGlobalValue (const GlobalOrdinal globalRow, const Scalar &value) |
Replace current value at the specified location with specified value. More... | |
void | sumIntoGlobalValue (const GlobalOrdinal globalRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault) |
Add value to existing value, using global (row) index. More... | |
void | replaceLocalValue (const LocalOrdinal myRow, const Scalar &value) |
Replace current value at the specified location with specified values. More... | |
void | sumIntoLocalValue (const LocalOrdinal myRow, const Scalar &value, const bool atomic=base_type::useAtomicUpdatesByDefault) |
Add value to existing value, using local (row) index. More... | |
Extraction methods | |
void | get1dCopy (const Teuchos::ArrayView< Scalar > &A) const |
Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management classes). More... | |
Teuchos::ArrayRCP< Scalar > | getDataNonConst () |
View of the local values of this vector. More... | |
Teuchos::ArrayRCP< const Scalar > | getData () const |
Const view of the local values of this vector. More... | |
Teuchos::RCP< const Vector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | offsetView (const Teuchos::RCP< const map_type > &subMap, const size_t offset) const |
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | offsetViewNonConst (const Teuchos::RCP< const map_type > &subMap, const size_t offset) |
Mathematical methods | |
dot_type | dot (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &y) const |
Return the dot product of this Vector and the input Vector x. More... | |
mag_type | norm1 () const |
Return the one-norm of this Vector. More... | |
mag_type | norm2 () const |
Return the two-norm of this Vector. More... | |
mag_type | normInf () const |
Return the infinity-norm of this Vector. More... | |
Scalar | meanValue () const |
Compute mean (average) value of this Vector. More... | |
Implementation of the Teuchos::Describable interface | |
virtual std::string | description () const |
Return a one-line description of this object. More... | |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Describe this object in a human-readable way to the given output stream. More... | |
Constructors and destructor | |
void | swap (MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &mv) |
Swap contents of mv with contents of *this . More... | |
Get a copy or view of a subset of rows and/or columns | |
The following methods get either a (deep) copy or a view (shallow copy) of a subset of rows and/or columns of the MultiVector. They return one of the following:
We prefer use of Kokkos classes to Teuchos Memory Management Classes. In particular, Teuchos::ArrayRCP reference counts are not thread safe, while Kokkos::View (and Kokkos::DualView) reference counts are thread safe. Not all of these methods are valid for a particular MultiVector. For instance, calling a method that accesses a view of the data in a 1-D format (i.e., get1dView) requires that the target MultiVector have constant stride. This category of methods also includes sync(), modify(), and getLocalView(), which help MultiVector implement DualView semantics. | |
Teuchos::RCP< MultiVector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subCopy (const Teuchos::Range1D &colRng) const |
Return a MultiVector with copies of selected columns. More... | |
Teuchos::RCP< MultiVector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subCopy (const Teuchos::ArrayView< const size_t > &cols) const |
Return a MultiVector with copies of selected columns. More... | |
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subView (const Teuchos::Range1D &colRng) const |
Return a const MultiVector with const views of selected columns. More... | |
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subView (const Teuchos::ArrayView< const size_t > &cols) const |
Return a const MultiVector with const views of selected columns. More... | |
Teuchos::RCP< MultiVector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subViewNonConst (const Teuchos::Range1D &colRng) |
Return a MultiVector with views of selected columns. More... | |
Teuchos::RCP< MultiVector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | subViewNonConst (const Teuchos::ArrayView< const size_t > &cols) |
Return a MultiVector with views of selected columns. More... | |
Teuchos::RCP< const Vector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | getVector (const size_t j) const |
Return a Vector which is a const view of column j. More... | |
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | getVectorNonConst (const size_t j) |
Return a Vector which is a nonconst view of column j. More... | |
Teuchos::ArrayRCP< const Scalar > | getData (size_t j) const |
Const view of the local values in a particular vector of this multivector. More... | |
Teuchos::ArrayRCP< Scalar > | getDataNonConst (size_t j) |
View of the local values in a particular vector of this multivector. More... | |
void | get1dCopy (const Teuchos::ArrayView< Scalar > &A, const size_t LDA) const |
Fill the given array with a copy of this multivector's local values. More... | |
void | get2dCopy (const Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > &ArrayOfPtrs) const |
Fill the given array with a copy of this multivector's local values. More... | |
Teuchos::ArrayRCP< const Scalar > | get1dView () const |
Const persisting (1-D) view of this multivector's local values. More... | |
Teuchos::ArrayRCP < Teuchos::ArrayRCP< const Scalar > > | get2dView () const |
Return const persisting pointers to values. More... | |
Teuchos::ArrayRCP< Scalar > | get1dViewNonConst () |
Nonconst persisting (1-D) view of this multivector's local values. More... | |
Teuchos::ArrayRCP < Teuchos::ArrayRCP< Scalar > > | get2dViewNonConst () |
Return non-const persisting pointers to values. More... | |
dual_view_type::t_host::const_type | getLocalViewHost (Access::ReadOnlyStruct) const |
Return a read-only, up-to-date view of this MultiVector's local data on host. This requires that there are no live device-space views. More... | |
dual_view_type::t_host | getLocalViewHost (Access::ReadWriteStruct) |
Return a mutable, up-to-date view of this MultiVector's local data on host. This requires that there are no live device-space views. More... | |
dual_view_type::t_host | getLocalViewHost (Access::OverwriteAllStruct) |
Return a mutable view of this MultiVector's local data on host, assuming all existing data will be overwritten. This requires that there are no live device-space views. More... | |
dual_view_type::t_dev::const_type | getLocalViewDevice (Access::ReadOnlyStruct) const |
Return a read-only, up-to-date view of this MultiVector's local data on device. This requires that there are no live host-space views. More... | |
dual_view_type::t_dev | getLocalViewDevice (Access::ReadWriteStruct) |
Return a mutable, up-to-date view of this MultiVector's local data on device. This requires that there are no live host-space views. More... | |
dual_view_type::t_dev | getLocalViewDevice (Access::OverwriteAllStruct) |
Return a mutable view of this MultiVector's local data on device, assuming all existing data will be overwritten. This requires that there are no live host-space views. More... | |
wrapped_dual_view_type | getWrappedDualView () const |
Return the wrapped dual view holding this MultiVector's local data. More... | |
bool | need_sync () const |
Whether this MultiVector needs synchronization to the given space. More... | |
bool | need_sync_host () const |
Whether this MultiVector needs synchronization to the host. More... | |
bool | need_sync_device () const |
Whether this MultiVector needs synchronization to the device. More... | |
std::remove_reference < decltype(std::declval < dual_view_type >).template view< TargetDeviceType > ))>::type::const_type | getLocalView (Access::ReadOnlyStruct s) const |
Return a view of the local data on a specific device, with the given access mode. The return type is either dual_view_type::t_dev, dual_view_type::t_host, or the const_type of one of those. More... | |
std::remove_reference < decltype(std::declval < dual_view_type >).template view< TargetDeviceType > ))>::type | getLocalView (Access::ReadWriteStruct s) |
std::remove_reference < decltype(std::declval < dual_view_type >).template view< TargetDeviceType > ))>::type | getLocalView (Access::OverwriteAllStruct s) |
Mathematical methods | |
void | dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< dot_type > &dots) const |
Compute the dot product of each corresponding pair of vectors (columns) in A and B. More... | |
std::enable_if< !(std::is_same < dot_type, T >::value), void > ::type | dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< T > &dots) const |
Compute the dot product of each corresponding pair of vectors (columns) in A and B. More... | |
std::enable_if< !(std::is_same < dot_type, T >::value), void > ::type | dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, std::vector< T > &dots) const |
Like the above dot() overload, but for std::vector output. More... | |
void | dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Kokkos::View< dot_type *, Kokkos::HostSpace > &norms) const |
Compute the dot product of each corresponding pair of vectors (columns) in A and B, storing the result in a device View. More... | |
void | dot (typename std::enable_if< std::is_same< typename ViewType::value_type, dot_type >::value &&std::is_same< typename ViewType::memory_space, typename device_type::memory_space >::value, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >>::type &A, const ViewType &dots) const |
std::enable_if< !(std::is_same < dot_type, T >::value), void > ::type | dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Kokkos::View< T *, device_type > &dots) const |
Compute the dot product of each corresponding pair of vectors (columns) in A and B, storing the result in a device view. More... | |
void | abs (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) |
Put element-wise absolute values of input Multi-vector in target: A = abs(this) More... | |
void | reciprocal (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) |
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j). More... | |
void | scale (const Scalar &alpha) |
Scale in place: this = alpha*this . More... | |
void | scale (const Teuchos::ArrayView< const Scalar > &alpha) |
Scale each column in place: this[j] = alpha[j]*this[j] . More... | |
void | scale (const Kokkos::View< const impl_scalar_type *, device_type > &alpha) |
Scale each column in place: this[j] = alpha[j]*this[j] . More... | |
void | scale (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) |
Scale in place: this = alpha * A . More... | |
void | update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) |
Update: this = beta*this + alpha*A . More... | |
void | update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma) |
Update: this = gamma*this + alpha*A + beta*B . More... | |
void | norm1 (const Kokkos::View< mag_type *, Kokkos::HostSpace > &norms) const |
Compute the one-norm of each vector (column), storing the result in a host view. More... | |
std::enable_if< std::is_same < typename ViewType::value_type, mag_type > ::value &&std::is_same < typename ViewType::memory_space, typename device_type::memory_space > ::value >::type | norm1 (const ViewType &norms) const |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | norm1 (const Kokkos::View< T *, device_type > &norms) const |
Compute the one-norm of each vector (column), storing the result in a device view. More... | |
void | norm1 (const Teuchos::ArrayView< mag_type > &norms) const |
Compute the one-norm of each vector (column). More... | |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | norm1 (const Teuchos::ArrayView< T > &norms) const |
Compute the one-norm of each vector (column). More... | |
void | norm2 (const Kokkos::View< mag_type *, Kokkos::HostSpace > &norms) const |
Compute the two-norm of each vector (column), storing the result in a host View. More... | |
std::enable_if< std::is_same < typename ViewType::value_type, mag_type > ::value &&std::is_same < typename ViewType::memory_space, typename device_type::memory_space > ::value >::type | norm2 (const ViewType &norms) const |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | norm2 (const Kokkos::View< T *, device_type > &norms) const |
Compute the two-norm of each vector (column), storing the result in a device view. More... | |
void | norm2 (const Teuchos::ArrayView< mag_type > &norms) const |
Compute the two-norm of each vector (column). More... | |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | norm2 (const Teuchos::ArrayView< T > &norms) const |
Compute the two-norm of each vector (column). More... | |
void | normInf (const Kokkos::View< mag_type *, Kokkos::HostSpace > &norms) const |
Compute the infinity-norm of each vector (column), storing the result in a host View. More... | |
std::enable_if< std::is_same < typename ViewType::value_type, mag_type > ::value &&std::is_same < typename ViewType::memory_space, typename device_type::memory_space > ::value >::type | normInf (const ViewType &norms) const |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | normInf (const Kokkos::View< T *, device_type > &norms) const |
Compute the infinity-norm of each vector (column), storing the result in a device view. More... | |
void | normInf (const Teuchos::ArrayView< mag_type > &norms) const |
Compute the infinity-norm of each vector (column), storing the result in a Teuchos::ArrayView. More... | |
std::enable_if< !(std::is_same < mag_type, T >::value), void > ::type | normInf (const Teuchos::ArrayView< T > &norms) const |
Compute the infinity-norm of each vector (column), storing the result in a Teuchos::ArrayView. More... | |
void | meanValue (const Teuchos::ArrayView< impl_scalar_type > &means) const |
Compute mean (average) value of each column. More... | |
std::enable_if<!std::is_same < impl_scalar_type, T >::value, void >::type | meanValue (const Teuchos::ArrayView< T > &means) const |
void | multiply (Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta) |
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B) . More... | |
void | elementWiseMultiply (Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis) |
Multiply a Vector A elementwise by a MultiVector B. More... | |
Attribute access functions | |
size_t | getNumVectors () const |
Number of columns in the multivector. More... | |
size_t | getLocalLength () const |
Local number of rows on the calling process. More... | |
global_size_t | getGlobalLength () const |
Global number of rows in the multivector. More... | |
size_t | getStride () const |
Stride between columns in the multivector. More... | |
bool | isConstantStride () const |
Whether this multivector has constant stride between columns. More... | |
bool | aliases (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &other) const |
Whether this multivector's memory might alias other. This is conservative: if either this or other is not constant stride, then it simply checks whether the contiguous memory allocations overlap. It doesn't check whether the sets of columns overlap. This is a symmetric relation: X.aliases(Y) == Y.aliases(X). More... | |
Methods for use only by experts | |
virtual void | removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) |
Remove processes which contain no entries in this object's Map. More... | |
Public methods for redistributing data | |
void | doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
Import data into this object using an Import object ("forward mode"). More... | |
void | doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
Import data into this object using an Export object ("reverse mode"). More... | |
void | doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
Export data into this object using an Export object ("forward mode"). More... | |
void | doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
Export data into this object using an Import object ("reverse mode"). More... | |
void | beginImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
void | beginImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
void | beginExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
void | beginExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
void | endImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
void | endImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
void | endExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false) |
void | endExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false) |
bool | transferArrived () const |
Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step. More... | |
Attribute accessor methods | |
bool | isDistributed () const |
Whether this is a globally distributed object. More... | |
virtual Teuchos::RCP< const map_type > | getMap () const |
The Map describing the parallel distribution of this object. More... | |
I/O methods | |
void | print (std::ostream &os) const |
Print this object to the given output stream. More... | |
Protected Member Functions | |
virtual void | doTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode) |
Redistribute data across (MPI) processes. More... | |
virtual bool | reallocArraysForNumPacketsPerLid (const size_t numExportLIDs, const size_t numImportLIDs) |
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary. More... | |
void | beginTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode) |
Implementation detail of doTransfer. More... | |
Misc. implementation details | |
std::string | descriptionImpl (const std::string &className) const |
Implementation of description() for this class, and its subclass Vector. More... | |
std::string | localDescribeToString (const Teuchos::EVerbosityLevel vl) const |
Print the calling process' verbose describe() information to the returned string. More... | |
void | describeImpl (Teuchos::FancyOStream &out, const std::string &className, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Implementation of describe() for this class, and its subclass Vector. More... | |
bool | vectorIndexOutOfRange (const size_t VectorIndex) const |
Teuchos::ArrayRCP< T > | getSubArrayRCP (Teuchos::ArrayRCP< T > arr, size_t j) const |
Persisting view of j-th column in the given ArrayRCP. More... | |
size_t | getOrigNumLocalRows () const |
"Original" number of rows in the (local) data. More... | |
size_t | getOrigNumLocalCols () const |
"Original" number of columns in the (local) data. More... | |
Protected Attributes | |
wrapped_dual_view_type | view_ |
The Kokkos::DualView containing the MultiVector's data. More... | |
Teuchos::Array< size_t > | whichVectors_ |
Indices of columns this multivector is viewing. More... | |
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node > | |
Teuchos::RCP< MultiVector < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | createMultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t numVectors) |
Nonmember MultiVector "constructor": Create a MultiVector from a given Map. More... | |
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node > | |
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > | createCopy (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &src) |
Return a deep copy of the given Vector. More... | |
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node > | |
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | createVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map) |
Nonmember Vector "constructor": Create a Vector from a given Map. More... | |
Methods implemented by subclasses and used by doTransfer(). | |
The doTransfer() method uses the subclass' implementations of these methods to implement data transfer. Subclasses of DistObject must implement these methods. This is an instance of the Template Method Pattern. ("Template" here doesn't mean "C++ template"; it means "pattern with holes that are filled in by the subclass' method implementations.") | |
Teuchos::RCP< const map_type > | map_ |
The Map over which this object is distributed. More... | |
Kokkos::DualView< packet_type *, buffer_device_type > | imports_ |
Buffer into which packed data are imported (received from other processes). More... | |
Kokkos::DualView< size_t *, buffer_device_type > | numImportPacketsPerLID_ |
Number of packets to receive for each receive operation. More... | |
Kokkos::DualView< packet_type *, buffer_device_type > | exports_ |
Buffer from which packed data are exported (sent to other processes). More... | |
Kokkos::DualView< size_t *, buffer_device_type > | numExportPacketsPerLID_ |
Number of packets to send for each send operation. More... | |
virtual void | copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) |
Perform copies and permutations that are local to the calling (MPI) process. More... | |
virtual void | copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM, const execution_space &space) |
Same as copyAndPermute, but do operations in space . More... | |
virtual void | packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) |
Pack data and metadata for communication (sends). More... | |
virtual void | packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets, const execution_space &space) |
Same as packAndPrepare, but in an execution space instance. More... | |
virtual void | unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode) |
Perform any unpacking and combining after communication. More... | |
virtual void | unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode, const execution_space &space) |
std::unique_ptr< std::string > | createPrefix (const char className[], const char methodName[]) const |
Post-construction modification routines | |
void | replaceGlobalValue (const GlobalOrdinal gblRow, const size_t col, const impl_scalar_type &value) |
Replace value in host memory, using global row index. More... | |
std::enable_if<!std::is_same < T, impl_scalar_type >::value &&std::is_convertible< T, impl_scalar_type >::value, void >::type | replaceGlobalValue (GlobalOrdinal globalRow, size_t col, const T &value) |
Like the above replaceGlobalValue, but only enabled if T differs from impl_scalar_type. More... | |
void | sumIntoGlobalValue (const GlobalOrdinal gblRow, const size_t col, const impl_scalar_type &value, const bool atomic=useAtomicUpdatesByDefault) |
Update (+=) a value in host memory, using global row index. More... | |
std::enable_if<!std::is_same < T, impl_scalar_type >::value &&std::is_convertible< T, impl_scalar_type >::value, void >::type | sumIntoGlobalValue (const GlobalOrdinal gblRow, const size_t col, const T &val, const bool atomic=useAtomicUpdatesByDefault) |
Like the above sumIntoGlobalValue, but only enabled if T differs from impl_scalar_type. More... | |
void | replaceLocalValue (const LocalOrdinal lclRow, const size_t col, const impl_scalar_type &value) |
Replace value in host memory, using local (row) index. More... | |
std::enable_if<!std::is_same < T, impl_scalar_type >::value &&std::is_convertible< T, impl_scalar_type >::value, void >::type | replaceLocalValue (const LocalOrdinal lclRow, const size_t col, const T &val) |
Like the above replaceLocalValue, but only enabled if T differs from impl_scalar_type. More... | |
void | sumIntoLocalValue (const LocalOrdinal lclRow, const size_t col, const impl_scalar_type &val, const bool atomic=useAtomicUpdatesByDefault) |
Update (+=) a value in host memory, using local row index. More... | |
std::enable_if<!std::is_same < T, impl_scalar_type >::value &&std::is_convertible< T, impl_scalar_type >::value, void >::type | sumIntoLocalValue (const LocalOrdinal lclRow, const size_t col, const T &val, const bool atomic=useAtomicUpdatesByDefault) |
Like the above sumIntoLocalValue, but only enabled if T differs from impl_scalar_type. More... | |
void | putScalar (const Scalar &value) |
Set all values in the multivector with the given value. More... | |
std::enable_if<!std::is_same < T, impl_scalar_type >::value &&std::is_convertible< T, impl_scalar_type >::value, void >::type | putScalar (const T &value) |
Set all values in the multivector with the given value. More... | |
void | randomize () |
Set all values in the multivector to pseudorandom numbers. More... | |
void | randomize (const Scalar &minVal, const Scalar &maxVal) |
Set all values in the multivector to pseudorandom numbers in the given range. More... | |
void | replaceMap (const Teuchos::RCP< const map_type > &map) |
Replace the underlying Map in place. More... | |
void | reduce () |
Sum values of a locally replicated multivector across all processes. More... | |
static const bool | useAtomicUpdatesByDefault |
Whether sumIntoLocalValue and sumIntoGlobalValue should use atomic updates by default. More... | |
Implementation of Tpetra::DistObject | |
bool | importsAreAliased () |
using | buffer_device_type = typename DistObject< scalar_type, local_ordinal_type, global_ordinal_type, node_type >::buffer_device_type |
Kokkos::DualView < impl_scalar_type *, buffer_device_type > | unaliased_imports_ |
virtual bool | checkSizes (const SrcDistObject &sourceObj) override |
Whether data redistribution between sourceObj and this object is legal. More... | |
virtual size_t | constantNumberOfPackets () const override |
Number of packets to send per LID. More... | |
virtual void | copyAndPermute (const SrcDistObject &sourceObj, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM, const execution_space &space) override |
virtual void | copyAndPermute (const SrcDistObject &sourceObj, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override |
virtual void | packAndPrepare (const SrcDistObject &sourceObj, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< impl_scalar_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type >, size_t &constantNumPackets, const execution_space &space) override |
virtual void | packAndPrepare (const SrcDistObject &sourceObj, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< impl_scalar_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type >, size_t &constantNumPackets) override |
virtual void | unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< impl_scalar_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type >, const size_t constantNumPackets, const CombineMode CM, const execution_space &space) override |
virtual void | unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< impl_scalar_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type >, const size_t constantNumPackets, const CombineMode CM) override |
virtual bool | reallocImportsIfNeeded (const size_t newSize, const bool verbose, const std::string *prefix, const bool areRemoteLIDsContiguous=false, const CombineMode CM=INSERT) override |
Reallocate imports_ if needed. More... | |
A distributed dense vector.
Scalar | The type of each entry of the vector. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double .) |
LocalOrdinal | The type of local indices. See the documentation of Map for requirements. |
GlobalOrdinal | The type of global indices. See the documentation of Map for requirements. |
Node | The Kokkos Node type. See the documentation of Map for requirements. |
This class inherits from MultiVector, and has the same template parameters. A Vector is a special case of a MultiVector that has only one vector (column). It may be used wherever a MultiVector may be used. Please see the documentation of MultiVector for more details.
Definition at line 44 of file Tpetra_Vector_decl.hpp.
typedef Scalar Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type |
This class' first template parameter; the type of each entry in the Vector.
Definition at line 55 of file Tpetra_Vector_decl.hpp.
typedef base_type::impl_scalar_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type |
The type used internally in place of Scalar
.
Some Scalar
types might not work with Kokkos on all execution spaces, due to missing CUDA device macros or volatile overloads. The C++ standard type std::complex<T> has this problem. To fix this, we replace std::complex<T> values internally with the (usually) bitwise identical type Kokkos::complex<T>. The latter is the impl_scalar_type
corresponding to Scalar
= std::complex.
Definition at line 65 of file Tpetra_Vector_decl.hpp.
typedef LocalOrdinal Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type |
This class' second template parameter; the type of local indices.
Definition at line 67 of file Tpetra_Vector_decl.hpp.
typedef GlobalOrdinal Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type |
This class' third template parameter; the type of global indices.
Definition at line 69 of file Tpetra_Vector_decl.hpp.
typedef Node::device_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::device_type |
The Kokkos device type.
Definition at line 71 of file Tpetra_Vector_decl.hpp.
typedef Node Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type |
The Kokkos Node type.
Definition at line 74 of file Tpetra_Vector_decl.hpp.
typedef base_type::dot_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot_type |
Type of an inner ("dot") product result.
This is usually the same as impl_scalar_type
, but may differ if impl_scalar_type
is e.g., an uncertainty quantification type from the Stokhos package.
Definition at line 81 of file Tpetra_Vector_decl.hpp.
typedef base_type::mag_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type |
Type of a norm result.
This is usually the same as the type of the magnitude (absolute value) of impl_scalar_type
, but may differ if impl_scalar_type
is e.g., an uncertainty quantification type from the Stokhos package.
Definition at line 89 of file Tpetra_Vector_decl.hpp.
typedef base_type::dual_view_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dual_view_type |
Kokkos::DualView specialization used by this class.
Definition at line 92 of file Tpetra_Vector_decl.hpp.
typedef base_type::wrapped_dual_view_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::wrapped_dual_view_type |
WrappedDualView specialization used by this class.
Definition at line 95 of file Tpetra_Vector_decl.hpp.
typedef base_type::map_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::map_type |
The type of the Map specialization used by this class.
Definition at line 98 of file Tpetra_Vector_decl.hpp.
|
inherited |
Type of the (new) Kokkos execution space.
The execution space implements parallel operations, like parallel_for, parallel_reduce, and parallel_scan.
Definition at line 411 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
The type of each datum being sent or received in an Import or Export.
Note that this type does not always correspond to the Scalar
template parameter of subclasses.
Definition at line 304 of file Tpetra_DistObject_decl.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | ) |
Default constructor: makes a Vector with no rows or columns.
Definition at line 26 of file Tpetra_Vector_def.hpp.
|
explicit |
Basic constructor.
map | [in] The Vector's Map. The Map describes the distribution of rows over process(es) in the Map's communicator. |
zeroOut | [in] If true (the default), require that all the Vector's entries be zero on return. If false, the Vector's entries have undefined values on return, and must be set explicitly. |
Definition at line 32 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | source, |
const Teuchos::DataAccess | copyOrView | ||
) |
Copy constructor (shallow or deep copy).
source | [in] The Vector to copy. |
copyOrView | [in] If Teuchos::View, return a shallow copy (a view) of source . If Teuchos::Copy, return a deep copy of source . Regardless, the result has "view semantics." This means that copy construction or assignment (operator=) with the resulting object will always do a shallow copy, and will transmit view semantics to the result of the shallow copy. |
Definition at line 39 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | source, |
const Teuchos::RCP< const map_type > & | map, | ||
const local_ordinal_type | rowOffset = 0 |
||
) |
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Teuchos::RCP< const map_type > & | map, |
const Teuchos::ArrayView< const Scalar > & | A | ||
) |
Set vector values from an existing array (copy)
Definition at line 54 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Teuchos::RCP< const map_type > & | map, |
const dual_view_type & | view | ||
) |
Expert mode constructor, that takes a Kokkos::DualView of the Vector's data, and returns a Vector that views those data.
See the documentation of the MultiVector (parent class) constructor that takes the same arguments.
map | [in] Map describing the distribution of rows. |
view | [in] View of the data (shallow copy). |
Definition at line 61 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Teuchos::RCP< const map_type > & | map, |
const dual_view_type & | view, | ||
const dual_view_type & | origView | ||
) |
Expert mode constructor, that takes a Kokkos::DualView of the Vector's data and the "original" Kokkos::DualView of the data, and returns a Vector that views those data.
See the documentation of the MultiVector (parent class) constructor that takes the same arguments.
map | [in] Map describing the distribution of rows. |
view | [in] View of the data (shallow copy). |
origView | [in] "Original" view of the data (shallow copy). |
Definition at line 68 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const Teuchos::RCP< const map_type > & | map, |
const wrapped_dual_view_type & | d_view | ||
) |
Expert mode constructor, that takes a WrappedDualView of the MultiVector's data.
Definition at line 77 of file Tpetra_Vector_def.hpp.
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector | ( | const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | X, |
const size_t | j | ||
) |
Create a Vector that views a single column of the input MultiVector.
X | [in] Input MultiVector to view (in possibly nonconst fashion). |
j | [in] The column of X to view. |
Definition at line 84 of file Tpetra_Vector_def.hpp.
|
default |
Copy constructor (shallow copy).
Vector's copy constructor always does a shallow copy. Use the nonmember function Tpetra::deep_copy
(see Tpetra_MultiVector_decl.hpp
) to deep-copy one existing Vector to another, and use the two-argument "copy constructor" below (with copyOrView=Teuchos::Copy
) to create a Vector that is a deep copy of an existing Vector.
|
default |
Move constructor (shallow move).
|
virtualdefault |
Destructor (virtual for memory safety of derived classes).
=default
declarations for copy construction, move construction, copy assignment, and move assignment.
|
default |
Copy assignment (shallow copy).
Vector's copy assignment operator always does a shallow copy. Use the nonmember function Tpetra::deep_copy
(see Tpetra_MultiVector_decl.hpp
) to deep-copy one existing Vector to another.
|
default |
Move assigment (shallow move).
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValue | ( | const GlobalOrdinal | globalRow, |
const Scalar & | value | ||
) |
Replace current value at the specified location with specified value.
globalRow
must be a valid global element on this node, according to the row map. Definition at line 92 of file Tpetra_Vector_def.hpp.
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValue | ( | const GlobalOrdinal | globalRow, |
const Scalar & | value, | ||
const bool | atomic = base_type::useAtomicUpdatesByDefault |
||
) |
Add value to existing value, using global (row) index.
Add the given value to the existing value at row globalRow
(a global index).
This method affects the host memory version of the data. If the DeviceType
template parameter is a device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.
globalRow | [in] Global row index of the entry to modify. This must be a valid global row index on the calling process with respect to the Vector's Map. |
value | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
Definition at line 99 of file Tpetra_Vector_def.hpp.
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValue | ( | const LocalOrdinal | myRow, |
const Scalar & | value | ||
) |
Replace current value at the specified location with specified values.
localRow
must be a valid local element on this node, according to the row map. Definition at line 109 of file Tpetra_Vector_def.hpp.
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValue | ( | const LocalOrdinal | myRow, |
const Scalar & | value, | ||
const bool | atomic = base_type::useAtomicUpdatesByDefault |
||
) |
Add value
to existing value, using local (row) index.
Add the given value to the existing value at row localRow
(a local index).
This method affects the host memory version of the data. If the DeviceType
template parameter is a device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.
localRow | [in] Local row index of the entry to modify. |
value | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
Definition at line 116 of file Tpetra_Vector_def.hpp.
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get1dCopy | ( | const Teuchos::ArrayView< Scalar > & | A | ) | const |
Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management classes).
Definition at line 126 of file Tpetra_Vector_def.hpp.
|
inline |
View of the local values of this vector.
Definition at line 306 of file Tpetra_Vector_decl.hpp.
|
inline |
Const view of the local values of this vector.
Definition at line 312 of file Tpetra_Vector_decl.hpp.
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot | ( | const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > & | y | ) | const |
Return the dot product of this Vector and the input Vector x.
Definition at line 134 of file Tpetra_Vector_def.hpp.
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm1 | ( | ) | const |
Return the one-norm of this Vector.
Definition at line 154 of file Tpetra_Vector_def.hpp.
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm2 | ( | ) | const |
Return the two-norm of this Vector.
Definition at line 164 of file Tpetra_Vector_def.hpp.
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normInf | ( | ) | const |
Return the infinity-norm of this Vector.
Definition at line 174 of file Tpetra_Vector_def.hpp.
Scalar Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::meanValue | ( | ) | const |
Compute mean (average) value of this Vector.
Definition at line 144 of file Tpetra_Vector_def.hpp.
|
virtual |
Return a one-line description of this object.
Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 184 of file Tpetra_Vector_def.hpp.
|
virtual |
Describe this object in a human-readable way to the given output stream.
You must call this method as a collective over all processes in this object's communicator.
out | [out] Output stream to which to write. Only Process 0 in this object's communicator may write to the output stream. |
verbLevel | [in] Verbosity level. This also controls whether this method does any communication. At verbosity levels higher (greater) than Teuchos::VERB_LOW, this method may behave as a collective over the object's communicator. |
Teuchos::FancyOStream wraps std::ostream. It adds features like tab levels. If you just want to wrap std::cout, try this:
Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 191 of file Tpetra_Vector_def.hpp.
|
inherited |
Swap contents of mv
with contents of *this
.
|
inherited |
Replace value in host memory, using global row index.
Replace the current value at row gblRow
(a global index) and column col
with the given value. The column index is zero based.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
This method does not have an "atomic" option like sumIntoGlobalValue. This is deliberate. Replacement is not commutative, unlike += (modulo rounding error). Concurrent calls to replaceGlobalValue on different threads that modify the same entry/ies have undefined results. (It's not just that one thread might win; it's that the value might get messed up.)
gblRow | [in] Global row index of the entry to modify. This must be a valid global row index on the calling process with respect to the MultiVector's Map. |
col | [in] Column index of the entry to modify. |
value | [in] Incoming value to add to the entry. |
|
inlineinherited |
Like the above replaceGlobalValue, but only enabled if T differs from impl_scalar_type.
This method only exists if its template parameter T
and impl_scalar_type differ, and if it is syntactically possible to convert T
to impl_scalar_type. This method is mainly useful for backwards compatibility, when the Scalar template parameter differs from impl_scalar_type. That is commonly only the case when Scalar is std::complex<U> for some type U.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
This method does not have an "atomic" option like sumIntoGlobalValue. This is deliberate. Replacement is not commutative, unlike += (modulo rounding error). Concurrent calls to replaceGlobalValue on different threads that modify the same entry/ies have undefined results. (It's not just that one thread might win; it's that the value might get messed up.)
gblRow | [in] Global row index of the entry to modify. This must be a valid global row index on the calling process with respect to the MultiVector's Map. |
col | [in] Column index of the entry to modify. |
value | [in] Incoming value to add to the entry. |
Definition at line 888 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Update (+=) a value in host memory, using global row index.
Add the given value to the existing value at row gblRow
(a global index) and column col
. The column index is zero based.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
gblRow | [in] Global row index of the entry to modify. This must be a valid global row index on the calling process with respect to the MultiVector's Map. |
col | [in] Column index of the entry to modify. |
value | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
|
inlineinherited |
Like the above sumIntoGlobalValue, but only enabled if T differs from impl_scalar_type.
This method only exists if its template parameter T
and impl_scalar_type differ, and if it is syntactically possible to convert T
to impl_scalar_type. This method is mainly useful for backwards compatibility, when the Scalar template parameter differs from impl_scalar_type. That is commonly only the case when Scalar is std::complex<U> for some type U.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
gblRow | [in] Global row index of the entry to modify. This must be a valid global row index on the calling process with respect to the MultiVector's Map. |
col | [in] Column index of the entry to modify. |
val | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
Definition at line 955 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Replace value in host memory, using local (row) index.
Replace the current value at row lclRow
(a local index) and column col
with the given value. The column index is zero based.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
This method does not have an "atomic" option like sumIntoLocalValue. This is deliberate. Replacement is not commutative, unlike += (modulo rounding error). Concurrent calls to replaceLocalValue on different threads that modify the same entry/ies have undefined results. (It's not just that one thread might win; it's that the value might get messed up.)
lclRow | [in] Local row index of the entry to modify. Must be a valid local index in this MultiVector's Map on the calling process. |
col | [in] Column index of the entry to modify. |
value | [in] Incoming value to add to the entry. |
|
inlineinherited |
Like the above replaceLocalValue, but only enabled if T differs from impl_scalar_type.
This method only exists if its template parameter T
and impl_scalar_type differ, and if it is syntactically possible to convert T
to impl_scalar_type. This method is mainly useful for backwards compatibility, when the Scalar template parameter differs from impl_scalar_type. That is commonly only the case when Scalar is std::complex<U> for some type U.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
This method does not have an "atomic" option like sumIntoLocalValue. This is deliberate. Replacement is not commutative, unlike += (modulo rounding error). Concurrent calls to replaceLocalValue on different threads that modify the same entry/ies have undefined results. (It's not just that one thread might win; it's that the value might get messed up.)
lclRow | [in] Local row index of the entry to modify. Must be a valid local index in this MultiVector's Map on the calling process. |
col | [in] Column index of the entry to modify. |
val | [in] Incoming value to add to the entry. |
Definition at line 1032 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Update (+=) a value in host memory, using local row index.
Add the given value to the existing value at row localRow
(a local index) and column col
. The column index is zero based.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
lclRow | [in] Local row index of the entry to modify. Must be a valid local index in this MultiVector's Map on the calling process. |
col | [in] Column index of the entry to modify. |
val | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
|
inlineinherited |
Like the above sumIntoLocalValue, but only enabled if T differs from impl_scalar_type.
This method only exists if its template parameter T
and impl_scalar_type differ, and if it is syntactically possible to convert T
to impl_scalar_type. This method is mainly useful for backwards compatibility, when the Scalar template parameter differs from impl_scalar_type. That is commonly only the case when Scalar is std::complex<U> for some type U.
This method affects the host memory version of the data. If device_type
is a Kokkos device that has two memory spaces, and you want to modify the non-host version of the data, you must access the device View directly by calling getLocalView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation. This method calls sync_host() before modifying host data, and modify_host() afterwards.
lclRow | [in] Local row index of the entry to modify. |
col | [in] Column index of the entry to modify. |
val | [in] Incoming value to add to the entry. |
atomic | [in] Whether to use an atomic update. If this class' execution space is not Kokkos::Serial, then this is true by default, else it is false by default. |
Definition at line 1097 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Set all values in the multivector with the given value.
|
inlineinherited |
Set all values in the multivector with the given value.
This method only exists if its template parameter T
and impl_scalar_type differ, and if it is syntactically possible to convert T
to impl_scalar_type. This method is mainly useful for backwards compatibility, when the Scalar template parameter differs from impl_scalar_type. That is commonly only the case when Scalar is std::complex<U> for some type U.
Definition at line 1118 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Set all values in the multivector to pseudorandom numbers.
|
inherited |
Set all values in the multivector to pseudorandom numbers in the given range.
|
inherited |
Replace the underlying Map in place.
map->isCompatible (this->getMap ())
. "Similar" means that the communicators have the same number of processes, though these need not be in the same order (have the same assignments of ranks) or represent the same communication contexts. It means the same thing as the MPI_SIMILAR return value of MPI_COMM_COMPARE. See MPI 3.0 Standard, Section 6.4.1.This method replaces this object's Map with the given Map. This relabels the rows of the multivector using the global IDs in the input Map. Thus, it implicitly applies a permutation, without actually moving data. If the new Map's communicator has more processes than the original Map's communicator, it "projects" the MultiVector onto the new Map by filling in missing rows with zeros. If the new Map's communicator has fewer processes than the original Map's communicator, the method "forgets about" any rows that do not exist in the new Map. (It mathematical terms, if one considers a MultiVector as a function from one vector space to another, this operation restricts the range.)
This method must always be called collectively on the communicator with the largest number of processes: either this object's current communicator (this->getMap()->getComm()
), or the new Map's communicator (map->getComm()
). If the new Map's communicator has fewer processes, then the new Map must be null on processes excluded from the original communicator, and the current Map must be nonnull on all processes. If the new Map has more processes, then it must be nonnull on all those processes, and the original Map must be null on those processes which are not in the new Map's communicator. (The latter case can only happen to a MultiVector to which a replaceMap() operation has happened before.)
this->getMap ()->getComm ()
). We reserve the right to do checking in debug mode that requires this method to be called collectively in order not to deadlock.
|
inherited |
Sum values of a locally replicated multivector across all processes.
|
inherited |
Return a MultiVector with copies of selected columns.
|
inherited |
Return a MultiVector with copies of selected columns.
|
inherited |
Return a const MultiVector with const views of selected columns.
|
inherited |
Return a const MultiVector with const views of selected columns.
|
inherited |
Return a MultiVector with views of selected columns.
|
inherited |
Return a MultiVector with views of selected columns.
|
inherited |
Return a Vector which is a const view of column j.
|
inherited |
Return a Vector which is a nonconst view of column j.
|
inherited |
Const view of the local values in a particular vector of this multivector.
|
inherited |
View of the local values in a particular vector of this multivector.
|
inherited |
Fill the given array with a copy of this multivector's local values.
A | [out] View of the array to fill. We consider A as a matrix with column-major storage. |
LDA | [in] Leading dimension of the matrix A. |
|
inherited |
Fill the given array with a copy of this multivector's local values.
ArrayOfPtrs | [out] Array of arrays, one for each column of the multivector. On output, we fill ArrayOfPtrs[j] with the data for column j of this multivector. |
|
inherited |
Const persisting (1-D) view of this multivector's local values.
This method assumes that the columns of the multivector are stored contiguously. If not, this method throws std::runtime_error.
|
inherited |
Return const persisting pointers to values.
|
inherited |
Nonconst persisting (1-D) view of this multivector's local values.
This method assumes that the columns of the multivector are stored contiguously. If not, this method throws std::runtime_error.
|
inherited |
Return non-const persisting pointers to values.
|
inherited |
Return a read-only, up-to-date view of this MultiVector's local data on host. This requires that there are no live device-space views.
|
inherited |
Return a mutable, up-to-date view of this MultiVector's local data on host. This requires that there are no live device-space views.
|
inherited |
Return a mutable view of this MultiVector's local data on host, assuming all existing data will be overwritten. This requires that there are no live device-space views.
|
inherited |
Return a read-only, up-to-date view of this MultiVector's local data on device. This requires that there are no live host-space views.
|
inherited |
Return a mutable, up-to-date view of this MultiVector's local data on device. This requires that there are no live host-space views.
|
inherited |
Return a mutable view of this MultiVector's local data on device, assuming all existing data will be overwritten. This requires that there are no live host-space views.
|
inherited |
Return the wrapped dual view holding this MultiVector's local data.
|
inlineinherited |
Whether this MultiVector needs synchronization to the given space.
Definition at line 1453 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Whether this MultiVector needs synchronization to the host.
|
inherited |
Whether this MultiVector needs synchronization to the device.
|
inlineinherited |
Return a view of the local data on a specific device, with the given access mode. The return type is either dual_view_type::t_dev, dual_view_type::t_host, or the const_type of one of those.
TargetDeviceType | The Kokkos Device type whose data to return. |
For example, suppose you create a Tpetra::MultiVector for the Kokkos::Cuda device, like this:
If you want to get the CUDA device Kokkos::View as read-write, do this:
and if you want to get the host mirror of that View, do this:
Definition at line 1493 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
The "dot product" is the standard Euclidean inner product. If the type of entries of the vectors (impl_scalar_type) is complex, then A is transposed, not *this
. For example, if x and y each have one column, then x.dot (y, dots)
computes .
*this
and A have the same number of columns (vectors). dots
has at least as many entries as the number of columns in A.dots[j] == (this->getVector[j])->dot (* (A.getVector[j]))
|
inlineinherited |
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
T | The output type of the dot products. |
This method only exists if dot_type and T are different types. For example, if impl_scalar_type and dot_type differ, then this method ensures backwards compatibility with the previous interface (that returned dot products as impl_scalar_type rather than as dot_type). The complicated enable_if
expression just ensures that the method only exists if dot_type and T are different types; the method still returns void
, as above.
Definition at line 1547 of file Tpetra_MultiVector_decl.hpp.
|
inlineinherited |
Like the above dot() overload, but for std::vector output.
Definition at line 1562 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the dot product of each corresponding pair of vectors (columns) in A and B, storing the result in a device View.
The "dot product" is the standard Euclidean inner product. If the type of entries of the vectors (impl_scalar_type) is complex, then A is transposed, not *this
. For example, if x and y each have one column, then x.dot (y, dots)
computes .
A | [in] MultiVector with which to dot *this . |
dots | [out] Device View with getNumVectors() entries. |
this->getNumVectors () == A.getNumVectors ()
dots.extent (0) == A.getNumVectors ()
dots(j) == (this->getVector[j])->dot (* (A.getVector[j]))
|
inlineinherited |
Compute the dot product of each corresponding pair of vectors (columns) in A and B, storing the result in a device view.
T | The output type of the dot products. |
This method only exists if dot_type and T are different types. For example, if Scalar and dot_type differ, then this method ensures backwards compatibility with the previous interface (that returned dot products as Scalar rather than as dot_type). The complicated enable_if
expression just ensures that the method only exists if dot_type and T are different types; the method still returns void
, as above.
Definition at line 1621 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Put element-wise absolute values of input Multi-vector in target: A = abs(this)
|
inherited |
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
|
inherited |
Scale in place: this = alpha*this
.
Replace this MultiVector with alpha times this MultiVector. This method will always multiply, even if alpha is zero. That means, for example, that if *this
contains NaN entries before calling this method, the NaN entries will remain after this method finishes.
|
inherited |
Scale each column in place: this[j] = alpha[j]*this[j]
.
Replace each column j of this MultiVector with alpha[j]
times the current column j of this MultiVector. This method will always multiply, even if all the entries of alpha are zero. That means, for example, that if *this
contains NaN entries before calling this method, the NaN entries will remain after this method finishes.
|
inherited |
Scale each column in place: this[j] = alpha[j]*this[j]
.
Replace each column j of this MultiVector with alpha[j]
times the current column j of this MultiVector. This method will always multiply, even if all the entries of alpha are zero. That means, for example, that if *this
contains NaN entries before calling this method, the NaN entries will remain after this method finishes.
|
inherited |
Scale in place: this = alpha * A
.
Replace this MultiVector with scaled values of A. This method will always multiply, even if alpha is zero. That means, for example, that if *this
contains NaN entries before calling this method, the NaN entries will remain after this method finishes. It is legal for the input A to alias this MultiVector.
|
inherited |
Update: this = beta*this + alpha*A
.
Update this MultiVector with scaled values of A. If beta is zero, overwrite *this
unconditionally, even if it contains NaN entries. It is legal for the input A to alias this MultiVector.
|
inherited |
Update: this = gamma*this + alpha*A + beta*B
.
Update this MultiVector with scaled values of A and B. If gamma is zero, overwrite *this
unconditionally, even if it contains NaN entries. It is legal for the inputs A or B to alias this MultiVector.
|
inherited |
Compute the one-norm of each vector (column), storing the result in a host view.
norms | [out] Host View with getNumVectors() entries. |
norms.extent (0) == this->getNumVectors ()
norms(j) == (this->getVector[j])->norm1 (* (A.getVector[j]))
The one-norm of a vector is the sum of the magnitudes of the vector's entries. On exit, norms(j) is the one-norm of column j of this MultiVector.
|
inlineinherited |
Compute the one-norm of each vector (column), storing the result in a device view.
T | The output type of the dot products. |
See the above norm1() method for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms products as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1755 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the one-norm of each vector (column).
See the uppermost norm1() method above for documentation.
|
inlineinherited |
Compute the one-norm of each vector (column).
T | The output type of the norms. |
See the uppermost norm1() method above for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1791 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the two-norm of each vector (column), storing the result in a host View.
norms | [out] Host View with getNumVectors() entries. |
norms.extent (0) == this->getNumVectors ()
norms(j) == (this->getVector[j])->dot (* (A.getVector[j]))
The two-norm of a vector is the standard Euclidean norm, the square root of the sum of squares of the magnitudes of the vector's entries. On exit, norms(k) is the two-norm of column k of this MultiVector.
|
inlineinherited |
Compute the two-norm of each vector (column), storing the result in a device view.
See the above norm2() method for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1847 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the two-norm of each vector (column).
See the uppermost norm2() method above for documentation.
|
inlineinherited |
Compute the two-norm of each vector (column).
T | The output type of the norms. |
See the uppermost norm2() method above for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms products as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1883 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the infinity-norm of each vector (column), storing the result in a host View.
The infinity-norm of a vector is the maximum of the magnitudes of the vector's entries. On exit, norms(j) is the infinity-norm of column j of this MultiVector.
|
inlineinherited |
Compute the infinity-norm of each vector (column), storing the result in a device view.
See the above normInf() method for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1932 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute the infinity-norm of each vector (column), storing the result in a Teuchos::ArrayView.
See the uppermost normInf() method above for documentation.
|
inlineinherited |
Compute the infinity-norm of each vector (column), storing the result in a Teuchos::ArrayView.
T | The output type of the norms. |
See the uppermost normInf() method above for documentation.
This method only exists if mag_type and T are different types. For example, if Teuchos::ScalarTraits<Scalar>::magnitudeType and mag_type differ, then this method ensures backwards compatibility with the previous interface (that returned norms products as Teuchos::ScalarTraits<Scalar>::magnitudeType rather than as mag_type). The complicated enable_if
expression just ensures that the method only exists if mag_type and T are different types; the method still returns void
, as above.
Definition at line 1970 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Compute mean (average) value of each column.
The outcome of this routine is undefined for non-floating point scalar types (e.g., int).
|
inherited |
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B)
.
If beta is zero, overwrite *this
unconditionally, even if it contains NaN entries. This imitates the semantics of analogous BLAS routines like DGEMM.
|
inherited |
Multiply a Vector A elementwise by a MultiVector B.
Compute this = scalarThis * this + scalarAB * B @ A
where </tt> denotes element-wise multiplication. In pseudocode, if C denotes
*this
MultiVector:
for all rows i and columns j of C.
B must have the same dimensions as
*this
, while A must have the same number of rows but a single column.
We do not require that A, B, and
*this
have compatible Maps, as long as the number of rows in A, B, and *this
on each process is the same. For example, one or more of these vectors might have a locally replicated Map, or a Map with a local communicator (MPI_COMM_SELF
). This case may occur in block relaxation algorithms when applying a diagonal scaling.
|
inherited |
Number of columns in the multivector.
|
inherited |
Local number of rows on the calling process.
|
inherited |
Global number of rows in the multivector.
|
inherited |
Stride between columns in the multivector.
This is only meaningful if isConstantStride()
returns true.
|
inherited |
Whether this multivector has constant stride between columns.
|
inherited |
Whether this multivector's memory might alias other. This is conservative: if either this or other is not constant stride, then it simply checks whether the contiguous memory allocations overlap. It doesn't check whether the sets of columns overlap. This is a symmetric relation: X.aliases(Y) == Y.aliases(X).
|
overridevirtualinherited |
Remove processes owning zero rows from the Map and their communicator.
newMap | [in] This must be the result of calling the removeEmptyProcesses() method on the row Map. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes. |
|
virtualinherited |
Remove processes which contain no entries in this object's Map.
On input, this object is distributed over the Map returned by getMap() (the "original Map," with its communicator, the "original communicator"). The input newMap
of this method must be the same as the result of calling getMap()->removeEmptyProcesses()
. On processes in the original communicator which contain zero entries ("excluded processes," as opposed to "included processes"), the input newMap
must be Teuchos::null
(which is what getMap()->removeEmptyProcesses()
returns anyway).
On included processes, reassign this object's Map (that would be returned by getMap()) to the input newMap
, and do any work that needs to be done to restore correct semantics. On excluded processes, free any data that needs freeing, and do any other work that needs to be done to restore correct semantics.
This method has collective semantics over the original communicator. On exit, the only method of this object which is safe to call on excluded processes is the destructor. This implies that subclasses' destructors must not contain communication operations.
|
inlineinherited |
Set whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.
Definition at line 2139 of file Tpetra_MultiVector_decl.hpp.
|
inlineinherited |
Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.
Definition at line 2156 of file Tpetra_MultiVector_decl.hpp.
|
inherited |
Copy the contents of src
into *this
(deep copy).
src | [in] Source MultiVector (input of the deep copy). |
! src.getMap ().is_null () && ! this->getMap ().is_null ()
src.getMap ()->isCompatible (* (this->getMap ())
src
or *this
remain valid.*this
, or otherwise change its dimensions. This is not an assignment operator; it does not change anything in *this
other than the contents of storage.
|
inherited |
Return another MultiVector with the same entries, but converted to a different Scalar type T
.
|
inherited |
src | [in] MultiVector |
! vec.getMap ().is_null () && ! this->getMap ().is_null ()
vec.getMap ()->isCompatible (* (this->getMap ())
src
or *this
remain valid.
|
protectedinherited |
Implementation of description() for this class, and its subclass Vector.
className | [in] Name of the class calling this method: Either "Tpetra::MultiVector" or "Tpetra::Vector" (no quotes in the string, in either case). |
|
protectedinherited |
Print the calling process' verbose describe() information to the returned string.
This is an implementation detail of describe().
vl | [in] Verbosity level with which to print. |
|
protectedinherited |
Implementation of describe() for this class, and its subclass Vector.
out | [out] Output stream to which to write. Only Process 0 in this object's communicator may write to the output stream. |
className | [in] Name of the class calling this method. |
verbLevel | [in] Verbosity level. This also controls whether this method does any communication. At verbosity levels higher (greater) than Teuchos::VERB_LOW, this method behaves as a collective over the object's communicator. |
|
protectedinherited |
Persisting view of j-th column in the given ArrayRCP.
This method considers isConstantStride(). The ArrayRCP may correspond either to a compute buffer or a host view.
|
protectedinherited |
"Original" number of rows in the (local) data.
|
protectedinherited |
"Original" number of columns in the (local) data.
|
overrideprotectedvirtualinherited |
Whether data redistribution between sourceObj
and this object is legal.
This method is called in DistObject::doTransfer() to check whether data redistribution between the two objects is legal.
Implements Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
overrideprotectedvirtualinherited |
Number of packets to send per LID.
Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
protectedvirtualinherited |
Perform copies and permutations that are local to the calling (MPI) process.
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, packs the source object's data.
permuteToLIDs.need_sync_host()
, permuteToLIDs.need_sync_device()
, permuteFromLIDs.need_sync_host()
, and permuteFromLIDs.need_sync_device()
are all false.source | [in] On entry, the source object of the Export or Import operation. |
numSameIDs | [in] The number of elements that are the same on the source and target objects. These elements live on the same process in both the source and target objects. |
permuteToLIDs | [in] List of the elements that are permuted. They are listed by their local index (LID) in the destination object. |
permuteFromLIDs | [in] List of the elements that are permuted. They are listed by their local index (LID) in the source object. |
CM | [in] CombineMode to be used during copyAndPermute; may or may not be used by the particular object being called; behavior with respect to CombineMode may differ by object. |
|
protectedvirtualinherited |
Same as copyAndPermute, but do operations in space
.
|
protectedvirtualinherited |
Pack data and metadata for communication (sends).
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, packs the source object's data.
exportLIDs.need_sync_host ()
and exportLIDs.need_sync_device()
are both false.source | [in] Source object for the redistribution. |
exportLIDs | [in] List of the entries (as local IDs in the source object) that Tpetra will send to other processes. |
exports | [out] On exit, the packed data to send. Implementations must reallocate this as needed (prefer reusing the existing allocation if possible), and may modify and/or sync this wherever they like. |
numPacketsPerLID | [out] On exit, the implementation of this method must do one of two things: either set numPacketsPerLID[i] to the number of packets to be packed for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation must not modify the entries of numPacketsPerLID . If the former, the implementation may sync numPacketsPerLID this wherever it likes, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
constantNumPackets | [out] On exit, 0 if the number of packets per LID could differ, else (if nonzero) the number of packets per LID (which must be constant). |
|
protectedvirtualinherited |
Same as packAndPrepare, but in an execution space instance.
|
protectedvirtualinherited |
Perform any unpacking and combining after communication.
Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this
, unpacks the received data into itself, possibly modifying its entries.
importLIDs.need_sync_host ()
and importLIDs.need_sync_device()
are both false.importLIDs | [in] List of the entries (as LIDs in the destination object) we received from other processes. |
imports | [in/out] On input: Buffer of received data to unpack. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
numPacketsPerLID | [in/out] On input: If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i]. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference. |
constantNumPackets | [in] If nonzero, then the number of packets per LID is the same for all entries ("constant") and constantNumPackets is that number. If zero, then numPacketsPerLID[i] is the number of packets to unpack for LID importLIDs[i] . |
combineMode | [in] The CombineMode to use when combining the imported entries with existing entries. |
|
overrideprotectedvirtualinherited |
Reallocate imports_ if needed.
This unfortunately must be declared protected, for the same reason that imports_ is declared protected.
newSize | [in] New size of imports_. |
verbose | [in] Whether to print verbose debugging output to stderr on every (MPI) process in the communicator. |
prefix | [in] If verbose is true , then this is a nonnull prefix to print at the beginning of each line of verbose debugging output. Otherwise, not used. |
We don't need a "reallocExportsIfNeeded" method, because exports_
always gets passed into packAndPrepare() by nonconst reference. Thus, that method can resize the DualView without needing to call other DistObject methods.
Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
inherited |
Import data into this object using an Import object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Import object if you want to do an Import, else use doExport() with a precomputed Export object.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
importer | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap() . |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Import data into this object using an Export object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doImport() that takes a precomputed Import object in that case.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
exporter | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap() . (Note the difference from forward mode.) |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Export data into this object using an Export object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Export object if you want to do an Export, else use doImport() with a precomputed Import object.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
exporter | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap() . |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Export data into this object using an Import object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this
object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doExport() that takes a precomputed Export object in that case.
"Restricted Mode" does two things:
*this
, in a "locallyFitted" sense. This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.
source | [in] The "source" object for redistribution. |
importer | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap() . (Note the difference from forward mode.) |
CM | [in] How to combine incoming data with the same global index. |
|
inherited |
Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step.
|
inherited |
Whether this is a globally distributed object.
For a definition of "globally distributed" (and its opposite, "locally replicated"), see the documentation of Map's isDistributed() method.
|
inlinevirtualinherited |
The Map describing the parallel distribution of this object.
Note that some Tpetra objects might be distributed using multiple Map objects. For example, CrsMatrix has both a row Map and a column Map. It is up to the subclass to decide which Map to use when invoking the DistObject constructor.
Definition at line 559 of file Tpetra_DistObject_decl.hpp.
|
inherited |
Print this object to the given output stream.
We generally assume that all MPI processes can print to the given stream.
|
protectedvirtualinherited |
Redistribute data across (MPI) processes.
src | [in] The source object, to redistribute into the target object, which is *this object. |
transfer | [in] The Export or Import object representing the communication pattern. (Details::Transfer is the common base class of these two objects.) |
modeString | [in] Human-readable string, for verbose debugging output and error output, explaining what function called this method. Example: "doImport (forward)", "doExport (reverse)". |
revOp | [in] Whether to do a forward or reverse mode redistribution. |
CM | [in] The combine mode that describes how to combine values that map to the same global ID on the same process. |
|
protectedvirtualinherited |
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.
numExportLIDs | [in] Number of entries in the exportLIDs input array argument of doTransfer(). |
numImportLIDs | [in] Number of entries in the remoteLIDs input array argument of doTransfer(). |
|
protectedinherited |
Implementation detail of doTransfer.
LID DualViews come from the Transfer object given to doTransfer. They are always sync'd on both host and device. Users must never attempt to modify or sync them.
|
related |
Nonmember MultiVector "constructor": Create a MultiVector from a given Map.
map | [in] Map describing the distribution of rows of the resulting MultiVector. |
numVectors | [in] Number of columns of the resulting MultiVector. |
|
related |
Return a deep copy of the given Vector.
|
staticprotectedinherited |
Whether sumIntoLocalValue and sumIntoGlobalValue should use atomic updates by default.
Definition at line 812 of file Tpetra_MultiVector_decl.hpp.
|
mutableprotectedinherited |
The Kokkos::DualView containing the MultiVector's data.
This has to be declared mutable
, so that get1dView() can retain its current const
marking, even though it has always implied a device->host synchronization. Lesson to the reader: Use const
sparingly!
Definition at line 2212 of file Tpetra_MultiVector_decl.hpp.
|
protectedinherited |
Indices of columns this multivector is viewing.
If this array has nonzero size, then this multivector is a view of another multivector (the "original" multivector). In that case, whichVectors_ contains the indices of the columns of the original multivector. Furthermore, isConstantStride() returns false in this case.
If this array has zero size, then this multivector is not a view of any other multivector. Furthermore, the stride between columns of this multivector is a constant: thus, isConstantStride() returns true.
Definition at line 2226 of file Tpetra_MultiVector_decl.hpp.
|
protectedinherited |
The Map over which this object is distributed.
Definition at line 968 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Buffer into which packed data are imported (received from other processes).
Unfortunately, I had to declare these protected, because CrsMatrix uses them at one point. Please, nobody else use them.
Definition at line 981 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Number of packets to receive for each receive operation.
This array is used in Distributor::doPosts() (and doReversePosts()) when starting the ireceive operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.
Definition at line 1021 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Buffer from which packed data are exported (sent to other processes).
Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.
Definition at line 1028 of file Tpetra_DistObject_decl.hpp.
|
protectedinherited |
Number of packets to send for each send operation.
This array is used in Distributor::doPosts() (and doReversePosts()) for preparing for the send operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Unfortunately, I had to declare this protected, because CrsMatrix uses them at one point. Please, nobody else use it.
Definition at line 1043 of file Tpetra_DistObject_decl.hpp.