Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Public Member Functions | Protected Attributes | List of all members
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false > Class Template Reference

One or more distributed dense vectors. More...

#include <Tpetra_KokkosRefactor_MultiVector_decl.hpp>

Inheritance diagram for Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >:
Inheritance graph
[legend]

Public Member Functions

virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &newMap)
 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_type > &src)
 Copy the contents of src into *this (deep copy). More...
 
Constructors and destructor
 MultiVector ()
 Default constructor: makes a MultiVector with no rows or columns. More...
 
 MultiVector (const Teuchos::RCP< const map_type > &map, const size_t numVecs, const bool zeroOut=true)
 Basic constuctor. More...
 
 MultiVector (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &source)
 Copy constructor (shallow copy!). More...
 
 MultiVector (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &source, const Teuchos::DataAccess copyOrView)
 Copy constructor, with option to do deep or shallow copy. More...
 
 MultiVector (const Teuchos::RCP< const map_type > &map, const Teuchos::ArrayView< const Scalar > &A, const size_t LDA, const size_t NumVectors)
 Create multivector by copying two-dimensional array of local data. More...
 
 MultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, const size_t NumVectors)
 Create multivector by copying array of views of local data. More...
 
 MultiVector (const Teuchos::RCP< const map_type > &map, const dual_view_type &view)
 Constructor, that takes a Kokkos::DualView of the MultiVector's data, and returns a MultiVector that views those data. More...
 
 MultiVector (const Teuchos::RCP< const map_type > &map, const typename dual_view_type::t_dev &d_view)
 Constructor, that takes a Kokkos::View of the MultiVector's data (living in the Device's memory space), and returns a MultiVector that views those data. More...
 
 MultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &map, const dual_view_type &view, const dual_view_type &origView)
 Expert mode constructor, that takes a Kokkos::DualView of the MultiVector's data and the "original" Kokkos::DualView of the data, and returns a MultiVector that views those data. More...
 
 MultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &map, const dual_view_type &view, const Teuchos::ArrayView< const size_t > &whichVectors)
 Expert mode constructor for noncontiguous views. More...
 
 MultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &map, const dual_view_type &view, const dual_view_type &origView, const Teuchos::ArrayView< const size_t > &whichVectors)
 Expert mode constructor for noncontiguous views, with original view. More...
 
template<class Node2 >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node2,
Node2::classic > > 
clone (const Teuchos::RCP< Node2 > &node2) const
 Return a deep copy of this MultiVector, with a different Node type. More...
 
virtual ~MultiVector ()
 Destructor (virtual for memory safety of derived classes). More...
 
Post-construction modification routines
void replaceGlobalValue (GlobalOrdinal globalRow, size_t col, const impl_scalar_type &value)
 Replace value, using global (row) index. More...
 
template<typename T >
Kokkos::Impl::enable_if
<!Kokkos::Impl::is_same< 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 (GlobalOrdinal globalRow, size_t col, const impl_scalar_type &value)
 Add value to existing value, using global (row) index. More...
 
template<typename T >
Kokkos::Impl::enable_if
<!Kokkos::Impl::is_same< T,
impl_scalar_type >::value,
void >::type 
sumIntoGlobalValue (GlobalOrdinal globalRow, size_t col, const T &value)
 Like the above sumIntoGlobalValue, but only enabled if T differs from impl_scalar_type. More...
 
void replaceLocalValue (LocalOrdinal localRow, size_t col, const impl_scalar_type &value)
 Replace value, using local (row) index. More...
 
template<typename T >
Kokkos::Impl::enable_if
<!Kokkos::Impl::is_same< T,
impl_scalar_type >::value,
void >::type 
replaceLocalValue (LocalOrdinal localRow, size_t col, const T &value)
 Like the above replaceLocalValue, but only enabled if T differs from impl_scalar_type. More...
 
void sumIntoLocalValue (LocalOrdinal localRow, size_t col, const impl_scalar_type &value)
 Add value to existing value, using local (row) index. More...
 
template<typename T >
Kokkos::Impl::enable_if
<!Kokkos::Impl::is_same< T,
impl_scalar_type >::value,
void >::type 
sumIntoLocalValue (LocalOrdinal localRow, size_t col, const T &value)
 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...
 
void randomize ()
 Set all values in the multivector to pseudorandom numbers. More...
 
void replaceMap (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &map)
 Replace the underlying Map in place. More...
 
void reduce ()
 Sum values of a locally replicated multivector across all processes. More...
 
MultiVector< Scalar,
LocalOrdinal, GlobalOrdinal,
node_type, false > & 
operator= (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type, false > &source)
 Assign the contents of source to this multivector (deep copy). More...
 
Data Copy and View get methods
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, node_type > > 
subCopy (const Teuchos::Range1D &colRng) const
 Return a MultiVector with copies of selected columns. More...
 
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, node_type > > 
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_type > > 
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_type > > 
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_type > > 
subViewNonConst (const Teuchos::Range1D &colRng)
 Return a MultiVector with views of selected columns. More...
 
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, node_type > > 
subViewNonConst (const Teuchos::ArrayView< const size_t > &cols)
 Return a MultiVector with views of selected columns. More...
 
Teuchos::RCP< const
MultiVector< Scalar,
LocalOrdinal, GlobalOrdinal,
node_type > > 
offsetView (const Teuchos::RCP< const map_type > &subMap, const size_t offset) const
 Return a const view of a subset of rows. More...
 
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, node_type > > 
offsetViewNonConst (const Teuchos::RCP< const map_type > &subMap, const size_t offset)
 Return a nonconst view of a subset of rows. More...
 
Teuchos::RCP< const Vector
< Scalar, LocalOrdinal,
GlobalOrdinal, node_type,
false > > 
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_type, false > > 
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...
 
KokkosClassic::MultiVector
< Scalar, node_type
getLocalMV () const
 A view of the underlying KokkosClassic::MultiVector object. More...
 
dual_view_type getDualView () const
 Get the Kokkos::DualView which implements local storage. More...
 
template<class TargetDeviceType >
void sync ()
 Update data on device or host only if data in the other space has been marked as modified. More...
 
template<class TargetDeviceType >
void modify ()
 Mark data as modified on the given device TargetDeviceType. More...
 
template<class TargetDeviceType >
Kokkos::Impl::if_c
< Kokkos::Impl::is_same
< typename
execution_space::memory_space,
typename
TargetDeviceType::memory_space >
::value, typename
dual_view_type::t_dev,
typename
dual_view_type::t_host >::type 
getLocalView () const
 Return a view of the local data on a specific device. More...
 
Mathematical methods
void dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A, const Teuchos::ArrayView< dot_type > &dots) const
 Compute the dot product of each corresponding pair of vectors (columns) in A and B. More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< dot_type, T >::value), void >
::type 
dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A, const Teuchos::ArrayView< T > &dots) const
 Compute the dot product of each corresponding pair of vectors (columns) in A and B. More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< dot_type, T >::value), void >
::type 
dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &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_type > &A, const Kokkos::View< dot_type *, execution_space > &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...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< dot_type, T >::value), void >
::type 
dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A, const Kokkos::View< T *, execution_space > &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_type > &A)
 Put element-wise absolute values of input Multi-vector in target: A = abs(this) More...
 
void reciprocal (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &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 *, execution_space > &alpha)
 Scale each column in place: this[j] = alpha[j]*this[j]. More...
 
void scale (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A)
 Scale in place: this = alpha * A. More...
 
void update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A, const Scalar &beta)
 Update: this = beta*this + alpha*A. More...
 
void update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &B, const Scalar &gamma)
 Update: this = gamma*this + alpha*A + beta*B. More...
 
void norm1 (const Kokkos::View< mag_type *, execution_space > &norms) const
 Compute the one-norm of each vector (column), storing the result in a device view. More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< mag_type, T >::value), void >
::type 
norm1 (const Kokkos::View< T *, execution_space > &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...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::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 *, execution_space > &norms) const
 Compute the two-norm of each vector (column), storing the result in a device view. More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< mag_type, T >::value), void >
::type 
norm2 (const Kokkos::View< T *, execution_space > &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...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::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 *, execution_space > &norms) const
 Compute the infinity-norm of each vector (column), storing the result in a device view. More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< mag_type, T >::value), void >
::type 
normInf (const Kokkos::View< T *, execution_space > &norms) const
 Compute the two-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). More...
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< mag_type, T >::value), void >
::type 
normInf (const Teuchos::ArrayView< T > &norms) const
 Compute the infinity-norm of each vector (column). More...
 
void normWeighted (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &weights, const Teuchos::ArrayView< mag_type > &norms) const
 
template<typename T >
Kokkos::Impl::enable_if
< !(Kokkos::Impl::is_same
< mag_type, T >::value), void >
::type 
normWeighted (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &weights, const Teuchos::ArrayView< T > &norms) const
 Compute the weighted 2-norm (RMS Norm) of each column. More...
 
void meanValue (const Teuchos::ArrayView< impl_scalar_type > &means) const
 Compute mean (average) value of each column. More...
 
template<typename T >
Kokkos::Impl::enable_if
<!Kokkos::Impl::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_type > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &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_type, false > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &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...
 
Overridden from Teuchos::Describable
virtual std::string description () const
 A simple one-line description of this object. More...
 
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print the object with the given verbosity level to a FancyOStream. More...
 

Protected Member Functions

Misc. implementation details
bool vectorIndexOutOfRange (const size_t VectorIndex) const
 
template<class T >
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...
 
Implementation of Tpetra::DistObject
virtual bool checkSizes (const SrcDistObject &sourceObj)
 Whether data redistribution between sourceObj and this object is legal. More...
 
virtual size_t constantNumberOfPackets () const
 Number of packets to send per LID. More...
 
virtual bool useNewInterface ()
 Whether this class implements the old or new interface of DistObject. More...
 
virtual void copyAndPermuteNew (const SrcDistObject &sourceObj, size_t numSameIDs, const Kokkos::View< const local_ordinal_type *, execution_space > &permuteToLIDs, const Kokkos::View< const local_ordinal_type *, execution_space > &permuteFromLIDs)
 
virtual void packAndPrepareNew (const SrcDistObject &sourceObj, const Kokkos::View< const local_ordinal_type *, execution_space > &exportLIDs, Kokkos::View< impl_scalar_type *, execution_space > &exports, const Kokkos::View< size_t *, execution_space > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor)
 
virtual void unpackAndCombineNew (const Kokkos::View< const LocalOrdinal *, execution_space > &importLIDs, const Kokkos::View< const impl_scalar_type *, execution_space > &imports, const Kokkos::View< size_t *, execution_space > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM)
 
void createViews () const
 
void createViewsNonConst (KokkosClassic::ReadWriteOption rwo)
 
void releaseViews () const
 

Protected Attributes

dual_view_type view_
 The Kokkos::DualView containing the MultiVector's data. More...
 
dual_view_type origView_
 The "original view" of the MultiVector's data. More...
 
Teuchos::Array< size_t > whichVectors_
 Indices of columns this multivector is viewing. More...
 

Typedefs to facilitate template metaprogramming.

typedef Scalar scalar_type
 This class' first template parameter; the type of each entry in the MultiVector. More...
 
typedef
Kokkos::Details::ArithTraits
< Scalar >::val_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
Kokkos::Compat::KokkosDeviceWrapperNode
< DeviceType > 
node_type
 This class' fourth template parameter; the Kokkos Node type. More...
 
typedef
Kokkos::Details::InnerProductSpaceTraits
< impl_scalar_type >::dot_type 
dot_type
 Type of an inner ("dot") product result. More...
 
typedef
Kokkos::Details::ArithTraits
< impl_scalar_type >::mag_type 
mag_type
 Type of a norm result. More...
 
typedef DeviceType execution_space
 Type of the (new) Kokkos execution space. More...
 
typedef Kokkos::DualView
< impl_scalar_type
**, Kokkos::LayoutLeft,
typename
execution_space::execution_space > 
dual_view_type
 Kokkos::DualView specialization used by this class. More...
 
typedef Map< LocalOrdinal,
GlobalOrdinal, node_type
map_type
 The type of the Map specialization used by this class. More...
 

Generic implementation of various norms

enum  EWhichNorm
 Input argument for normImpl() (which see). More...
 
void normImpl (const Kokkos::View< mag_type *, execution_space > &norms, const EWhichNorm whichNorm) const
 Compute the norm of each vector (column), storing the result in a device View. More...
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class DeviceType>
class Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >

One or more distributed dense vectors.

A "multivector" contains one or more dense vectors. All the vectors in a multivector have the same distribution of rows in parallel over the communicator used to create the multivector. Multivectors containing more than one vector are useful for algorithms that solve multiple linear systems at once, or that solve for a cluster of eigenvalues and their corresponding eigenvectors at once. These "block" algorithms often have accuracy or performance advantages over corresponding algorithms that solve for only one vector at a time. For example, working with multiple vectors at a time allows Tpetra to use faster BLAS 3 routines for local computations. It may also reduce the number of parallel reductions.

The Vector class implements the MultiVector interface, so if you only wish to work with a single vector at a time, you may simply use Vector instead of MultiVector. However, if you are writing solvers or preconditioners, you would do better to write to the MultiVector interface and always assume that each MultiVector contains more than one vector. This will make your solver or preconditioner more compatible with other Trilinos packages, and it will also let you exploit the performance optimizations mentioned above.

Template Parameters
ScalarThe type of each entry of the multivector. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double.)
LocalOrdinalThe type of local indices. See the documentation of Map for requirements.
GlobalOrdinalThe type of global indices. See the documentation of Map for requirements.
DeviceTypeThe Kokkos execution space type.

Prerequisites

Before reading the rest of this documentation, it helps to know a little bit about Kokkos. In particular, you should know about execution spaces, memory spaces, and shallow copy semantics. You should also know something about the Teuchos memory management classes, in particular Teuchos::RCP, though it helps to know a bit about Teuchos::ArrayRCP and Teuchos::ArrayView as well. You may also want to know about the differences between BLAS 1, 2, and 3 operations, and learn a little bit about MPI (the Message Passing Interface for distributed-memory programming). You won't have to use MPI directly to use MultiVector, but it helps to be familiar with the general idea of distributed storage of data over a communicator.

A MultiVector is a view of data

A MultiVector is a view of data. A view behaves like a pointer; it provides access to the original multivector's data without copying the data. This means that the copy constructor and assignment operator (operator=) do shallow copies. They do not copy the data; they just copy pointers and other "metadata." If you would like to copy a MultiVector into an existing MultiVector, call the nonmember function deep_copy(). If you would like to create a new MultiVector which is a deep copy of an existing MultiVector, call the nonmember function createCopy(), or use the two-argument copy constructor with Teuchos::Copy as the second argument.

Views have the additional property that they automatically handle deallocation. They use reference counting for this, much like how std::shared_ptr works. That means you do not have to worry about "freeing" a MultiVector after it has been created. Furthermore, you may pass shallow copies around without needing to worry about which is the "master" view of the data. There is no "master" view of the data; when the last view falls out of scope, the data will be deallocated.

This is what the documentation means when it speaks of view semantics. The opposite of that is copy or container semantics, where the copy constructor and operator= do deep copies (of the data). We say that "std::vector has container semantics," and "MultiVector has view semantics."

MultiVector also has "subview" methods that give results analogous to the Kokkos::subview() function. That is, they return a MultiVector which views some subset of another MultiVector's rows and columns. The subset of columns in a view need not be contiguous. For example, given a multivector X with 43 columns, it is possible to have a multivector Y which is a view of columns 1, 3, and 42 (zero-based indices) of X. We call such multivectors noncontiguous. They have the the property that isConstantStride() returns false.

Noncontiguous multivectors lose some performance advantages. For example, local computations may be slower, since Tpetra cannot use BLAS 3 routines (e.g., matrix-matrix multiply) on a noncontiguous multivectors without copying into temporary contiguous storage. For performance reasons, if you get a Kokkos::View of a noncontiguous MultiVector's local data, it does not correspond to the columns that the MultiVector views.

DualView semantics

Tpetra was designed to perform well on many different kinds of computers. Some computers have different memory spaces. For example, GPUs (Graphics Processing Units) by NVIDIA have "device memory" and "host memory." The GPU has faster access to device memory than host memory, but usually there is less device memory than host memory. Intel's "Knights Landing" architecture has two different memory spaces, also with different capacity and performance characteristics. Some architectures let the processor address memory in any space, possibly with a performance penalty. Others can only access data in certain spaces, and require a special system call to copy memory between spaces.

The Kokkos package provides abstractions for handling multiple memory spaces. In particular, Kokkos::DualView lets users "mirror" data that live in one space, with data in another space. It also lets users manually mark data in one space as modified (modify()), and synchronize (sync()) data from one space to another. The latter only actually copies if the data have been marked as modified. Users can access data in a particular space by calling view(). All three of these methods – modify(), sync(), and view() – are templated on the memory space. This is how users select the memory space on which they want the method to act.

MultiVector implements "DualView semantics." This means that it implements the above three operations:

If your computer only has one memory space, as with conventional single-core or multicore processors, you don't have to worry about this. You can ignore the modify() and sync() methods in that case.

How to access the local data

The getLocalView() method for getting a Kokkos::View, and getDualView() for getting a Kokkos::DualView, are the two main ways to access a MultiVector's local data. If you want to read or write the actual values in a multivector, this is what you want. The resulting Kokkos::View behaves like a 2-D array. You can address it using an index pair (i,j), where i is the local row index, and j is the column index.

MultiVector also has methods that return an Teuchos::ArrayRCP<Scalar> ("1-D view"), or a Teuchos::ArrayRCP<Teuchos::ArrayRCP<Scalar> > ("2-D view"). These exist only for backwards compatibility, and also give access to the local data.

All of these views only view local data. This means that the corresponding rows of the multivector are owned by the calling (MPI) process. You may not use these methods to access remote data, that is, rows that do not belong to the calling process.

MultiVector's public interface also has methods for modifying local data, like sumIntoLocalValue() and replaceGlobalValue(). These methods act on host data only. To access or modify device data, you must get the Kokkos::View and work with it directly.

Why won't you give me a raw pointer?

Tpetra was designed to allow different data representations underneath the same interface. This lets Tpetra run correctly and efficiently on many different kinds of hardware. These different kinds of hardware all have in common the following:

These conclusions have practical consequences for the MultiVector interface. In particular, we have deliberately made it difficult for you to access data directly by raw pointer. This is because the underlying layout may not be what you expect. The memory might not even be accessible from the host CPU. Instead, we give access through a Kokkos::View, which behaves like a 2-D array. You can ask the Kokkos::View for a raw pointer by calling its ptr_on_device() method, but then you are responsible for understanding its layout in memory.

Parallel distribution of data

A MultiVector's rows are distributed over processes in its (row) Map's communicator. A MultiVector is a DistObject; the Map of the DistObject tells which process in the communicator owns which rows. This means that you may use Import and Export operations to migrate between different distributions. Please refer to the documentation of Map, Import, and Export for more information.

MultiVector includes methods that perform parallel all-reduces. These include inner products and various kinds of norms. All of these methods have the same blocking semantics as MPI_Allreduce.

Warning
Some computational methods, such as inner products and norms, may return incorrect results if the MultiVector's Map is overlapping (not one-to-one) but not locally replicated. That is, if some but not all rows are shared by more than one process in the communicator, then inner products and norms may be wrong. This behavior may change in future releases.

Definition at line 286 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

Member Typedef Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Scalar Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::scalar_type

This class' first template parameter; the type of each entry in the MultiVector.

Definition at line 302 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Kokkos::Details::ArithTraits<Scalar>::val_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::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 missing volatile overloads of some methods. The C++ standard type std::complex<T> has this problem. To fix this, we replace std::complex<T> values internally with the bitwise identical type Kokkos::complex<T>. The latter is the impl_scalar_type corresponding to Scalar = std::complex<T>.

Most users don't need to know about this. Just be aware that if you ask for a Kokkos::View or Kokkos::DualView of the MultiVector's data, its entries have type impl_scalar_type, not scalar_type.

Definition at line 318 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef LocalOrdinal Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::local_ordinal_type

This class' second template parameter; the type of local indices.

Definition at line 320 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef GlobalOrdinal Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::global_ordinal_type

This class' third template parameter; the type of global indices.

Definition at line 322 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::node_type

This class' fourth template parameter; the Kokkos Node type.

Definition at line 324 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Kokkos::Details::InnerProductSpaceTraits<impl_scalar_type>::dot_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::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 336 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Kokkos::Details::ArithTraits<impl_scalar_type>::mag_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::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 344 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef DeviceType Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::execution_space

Type of the (new) Kokkos execution space.

The execution space implements parallel operations, like parallel_for, parallel_reduce, and parallel_scan. It also has a default memory space, in which the Tpetra object's data live.

Definition at line 352 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Kokkos::DualView<impl_scalar_type**, Kokkos::LayoutLeft, typename execution_space::execution_space> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dual_view_type

Kokkos::DualView specialization used by this class.

This is of interest to users who already have a Kokkos::DualView, and want the MultiVector to view it. By "view" it, we mean that the MultiVector doesn't copy the data in the DualView; it just hangs on to the pointer.

We take particular care to template the DualView on an execution space, rather than a memory space. This ensures that Tpetra will use exactly the specified execution space(s) and no others. This matters because View (and DualView) initialization is a parallel Kokkos kernel. If the View is templated on an execution space, Kokkos uses that execution space (and only that execution space) to initialize the View. This is what we want. If the View is templated on a memory space, Kokkos uses the memory space's default execution space to initialize. This is not necessarily what we want. For example, if building with OpenMP enabled, the default execution space for host memory is Kokkos::OpenMP, even if the user-specified DeviceType is Kokkos::Serial. That is why we go through the trouble of asking for the execution_space's execution space.

Definition at line 377 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
typedef Map<LocalOrdinal, GlobalOrdinal, node_type> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::map_type

The type of the Map specialization used by this class.

Definition at line 380 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

Member Enumeration Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
enum Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::EWhichNorm
protected

Input argument for normImpl() (which see).

Definition at line 1890 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( )

Default constructor: makes a MultiVector with no rows or columns.

Definition at line 252 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const map_type > &  map,
const size_t  numVecs,
const bool  zeroOut = true 
)

Basic constuctor.

Parameters
map[in] Map describing the distribution of rows.
numVecs[in] Number of vectors (columns).
zeroOut[in] Whether to initialize all the entries of the MultiVector to zero.

Definition at line 258 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  source)

Copy constructor (shallow copy!).

Definition at line 275 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  source,
const Teuchos::DataAccess  copyOrView 
)

Copy constructor, with option to do deep or shallow copy.

The Kokkos refactor version of Tpetra, unlike the "classic" version, always has view semantics. Thus, copyOrView = Teuchos::View has no effect, and copyOrView = Teuchos::Copy does not mark this MultiVector as having copy semantics. However, copyOrView = Teuchos::Copy will make the resulting MultiVector a deep copy of the input MultiVector.

Definition at line 286 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const map_type > &  map,
const Teuchos::ArrayView< const Scalar > &  A,
const size_t  LDA,
const size_t  NumVectors 
)

Create multivector by copying two-dimensional array of local data.

Parameters
map[in] The Map describing the distribution of rows of the multivector.
view[in] A view of column-major dense matrix data. The calling process will make a deep copy of this data.
LDA[in] The leading dimension (a.k.a. "stride") of the column-major input data.
NumVectors[in] The number of columns in the input data. This will be the number of vectors in the returned multivector.
Precondition
LDA >= A.size()
NumVectors > 0
Postcondition
isConstantStride() == true
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  map,
const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &  ArrayOfPtrs,
const size_t  NumVectors 
)

Create multivector by copying array of views of local data.

Parameters
map[in] The Map describing the distribution of rows of the multivector.
ArrayOfPtrs[in/out] Array of views of each column's data. The calling process will make a deep copy of this data.
NumVectors[in] The number of columns in the input data, and the number of elements in ArrayOfPtrs. This will be the number of vectors in the returned multivector.
Precondition
NumVectors > 0
NumVectors == ArrayOfPtrs.size()
Postcondition
constantStride() == true
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const map_type > &  map,
const dual_view_type view 
)

Constructor, that takes a Kokkos::DualView of the MultiVector's data, and returns a MultiVector that views those data.

To "view those data" means that this MultiVector and the input Kokkos::DualView point to the same data, just like two "raw" pointers (e.g., double*) can point to the same data. If you modify one, the other sees it (subject to the limitations of cache coherence).

Parameters
map[in] Map describing the distribution of rows.
view[in] Kokkos::DualView of the data to view.

Definition at line 318 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const map_type > &  map,
const typename dual_view_type::t_dev &  d_view 
)

Constructor, that takes a Kokkos::View of the MultiVector's data (living in the Device's memory space), and returns a MultiVector that views those data.

Parameters
map[in] Map describing the distribution of rows.
view[in] Kokkos::View of the data to view.

Q: What's the difference between this constructor (that takes a Kokkos::View), and the constructor above that takes a Kokkos::DualView?

A: Suppose that for the MultiVector's device type, there are actually two memory spaces (e.g., for Kokkos::Cuda with UVM off, assuming that this is allowed). In order for MultiVector to implement DualView semantics correctly, this constructor must allocate a Kokkos::View of host memory (or lazily allocate it on modify() or sync()).

Now suppose that you pass in the same Kokkos::View of device memory to two different MultiVector instances, X and Y. Each would allocate its own Kokkos::View of host memory. That means that X and Y would have different DualView instances, but their DualView instances would have the same device View.

Suppose that you do the following:

  1. Modify X on host (calling modify() correctly)
  2. Modify Y on host (calling modify() correctly)
  3. Sync Y to device (calling sync() correctly)
  4. Sync X to device (calling sync() correctly)

This would clobber Y's later changes on host, in favor of X's earlier changes on host. That could be very confusing. We allow this behavior because Kokkos::DualView allows it. That is, Kokkos::DualView also lets you get the device View, and hand it off to another Kokkos::DualView. It's confusing, but users need to know what they are doing if they start messing around with multiple memory spaces.

Definition at line 344 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  map,
const dual_view_type view,
const dual_view_type origView 
)

Expert mode constructor, that takes a Kokkos::DualView of the MultiVector's data and the "original" Kokkos::DualView of the data, and returns a MultiVector that views those data.

Warning
This constructor is only for expert users. We make no promises about backwards compatibility for this interface. It may change or go away at any time. It is mainly useful for Tpetra developers and we do not expect it to be useful for anyone else.
Parameters
map[in] Map describing the distribution of rows.
view[in] View of the data (shallow copy).
origView[in] The original view of the data.

The original view keeps the "original" dimensions. Doing so lets us safely construct a column Map view of a (domain Map view of a (column Map MultiVector)). The result of a Kokkos::subview does not remember the original dimensions of the view, and does not allow constructing a view with a superset of rows or columns, so we have to keep the original view.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  map,
const dual_view_type view,
const Teuchos::ArrayView< const size_t > &  whichVectors 
)

Expert mode constructor for noncontiguous views.

Warning
This constructor is only for expert users. We make no promises about backwards compatibility for this interface. It may change or go away at any time. It is mainly useful for Tpetra developers and we do not expect it to be useful for anyone else.

This constructor takes a Kokkos::DualView for the MultiVector to view, and a list of the columns to view, and returns a MultiVector that views those data. The resulting MultiVector does not have constant stride, that is, isConstantStride() returns false.

Parameters
map[in] Map describing the distribution of rows.
view[in] Device view to the data (shallow copy).
whichVectors[in] Which columns (vectors) to view.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::MultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  map,
const dual_view_type view,
const dual_view_type origView,
const Teuchos::ArrayView< const size_t > &  whichVectors 
)

Expert mode constructor for noncontiguous views, with original view.

Warning
This constructor is only for expert users. We make no promises about backwards compatibility for this interface. It may change or go away at any time. It is mainly useful for Tpetra developers and we do not expect it to be useful for anyone else.

This constructor takes a Kokkos::DualView for the MultiVector to view, a view of the "original" data, and a list of the columns to view, and returns a MultiVector that views those data. The resulting MultiVector does not have constant stride, that is, isConstantStride() returns false.

Parameters
map[in] Map describing the distribution of rows.
view[in] View of the data (shallow copy).
origView[in] The original view of the data.
whichVectors[in] Which columns (vectors) to view.

The original view keeps the "original" dimensions. Doing so lets us safely construct a column Map view of a (domain Map view of a (column Map MultiVector)). The result of a Kokkos::subview does not remember the original dimensions of the view, and does not allow constructing a view with a superset of rows or columns, so we have to keep the original view.

Definition at line 480 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::~MultiVector ( )
virtual

Destructor (virtual for memory safety of derived classes).

Definition at line 653 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<class Node2 >
Teuchos::RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2, Node2::classic> > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::clone ( const Teuchos::RCP< Node2 > &  node2) const

Return a deep copy of this MultiVector, with a different Node type.

Parameters
node2[in/out] The new Node type.
Warning
We prefer that you use Tpetra::deep_copy (see below) rather than this method. This method will go away at some point.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::replaceGlobalValue ( GlobalOrdinal  globalRow,
size_t  col,
const impl_scalar_type value 
)

Replace value, using global (row) index.

Replace the current value at row globalRow (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 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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Precondition
globalRow must be a valid global element on this process, according to the row Map.

Definition at line 3713 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if<! Kokkos::Impl::is_same<T, impl_scalar_type>::value, void>::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::replaceGlobalValue ( GlobalOrdinal  globalRow,
size_t  col,
const T &  value 
)
inline

Like the above replaceGlobalValue, but only enabled if T differs from impl_scalar_type.

This method only exists if the template parameter T and impl_scalar_type differ. If C++11 is enabled, we further require that it be possible to convert T to impl_scalar_type.

This method is mainly useful for backwards compatibility, when Scalar differs from impl_scalar_type.

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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Definition at line 648 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::sumIntoGlobalValue ( GlobalOrdinal  globalRow,
size_t  col,
const impl_scalar_type value 
)

Add value to existing value, using global (row) index.

Add the given value to the existing value at row globalRow (a global index) and column col. The column index is zero based.

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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Precondition
globalRow must be a valid global element on this process, according to the row Map.

Definition at line 3738 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if<! Kokkos::Impl::is_same<T, impl_scalar_type>::value, void>::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::sumIntoGlobalValue ( GlobalOrdinal  globalRow,
size_t  col,
const T &  value 
)
inline

Like the above sumIntoGlobalValue, but only enabled if T differs from impl_scalar_type.

This method only exists if the template parameter T and impl_scalar_type differ. If C++11 is enabled, we further require that it be possible to convert T to impl_scalar_type.

This method is mainly useful for backwards compatibility, when Scalar differs from impl_scalar_type.

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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Definition at line 699 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::replaceLocalValue ( LocalOrdinal  localRow,
size_t  col,
const impl_scalar_type value 
)

Replace value, using local (row) index.

Replace the current value at row localRow (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 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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Precondition
localRow must be a valid local element on this process, according to the row Map.

Definition at line 3651 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if<! Kokkos::Impl::is_same<T, impl_scalar_type>::value, void>::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::replaceLocalValue ( LocalOrdinal  localRow,
size_t  col,
const T &  value 
)
inline

Like the above replaceLocalValue, but only enabled if T differs from impl_scalar_type.

This method only exists if the template parameter T and impl_scalar_type differ. If C++11 is enabled, we further require that it be possible to convert T to impl_scalar_type.

This method is mainly useful for backwards compatibility, when Scalar differs from impl_scalar_type.

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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Definition at line 750 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::sumIntoLocalValue ( LocalOrdinal  localRow,
size_t  col,
const impl_scalar_type value 
)

Add value to existing value, 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 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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Precondition
localRow must be a valid local element on this process, according to the row Map.

Definition at line 3682 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if<! Kokkos::Impl::is_same<T, impl_scalar_type>::value, void>::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::sumIntoLocalValue ( LocalOrdinal  localRow,
size_t  col,
const T &  value 
)
inline

Like the above sumIntoLocalValue, but only enabled if T differs from impl_scalar_type.

This method only exists if the template parameter T and impl_scalar_type differ. If C++11 is enabled, we further require that it be possible to convert T to impl_scalar_type.

This method is mainly useful for backwards compatibility, when Scalar differs from impl_scalar_type.

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 DualView directly by calling getDualView(). Please see modify(), sync(), and the discussion of DualView semantics elsewhere in the documentation.

Definition at line 801 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::putScalar ( const Scalar &  value)

Set all values in the multivector with the given value.

Definition at line 1900 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::randomize ( )

Set all values in the multivector to pseudorandom numbers.

Note
The implementation of this method may depend on the Kokkos Node type and on what third-party libraries you have available. Do not expect repeatable results.
Behavior of this method may or may not depend on external use of the C library routines srand() and rand().
Warning
This method does not promise to use a distributed-memory parallel pseudorandom number generator. Corresponding values on different processes might be correlated. It also does not promise to use a high-quality pseudorandom number generator within each process.

Definition at line 1850 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::replaceMap ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  map)

Replace the underlying Map in place.

Warning
The normal use case of this method, with an input Map that is compatible with the object's current Map and has the same communicator, is safe. However, if the input Map has a different communicator (with a different number of processes, in particular) than this object's current Map, the semantics of this method are tricky. We recommend that only experts try the latter use case.
Precondition
If the new Map's communicator is similar to the original Map's communicator, then the original Map and new Map must be compatible: 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.
If the new Map's communicator contains more processes than the original Map's communicator, then the projection of the original Map onto the new communicator must be compatible with the new Map.
If the new Map's communicator contains fewer processes than the original Map's communicator, then the projection of the new Map onto the original communicator must be compatible with the original Map.

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.)

Warning
This method must always be called as a collective operation on all processes in the original communicator (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.
Note
This method does not do data redistribution. If you need to move data around, use Import or Export.

Definition at line 1983 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::reduce ( )

Sum values of a locally replicated multivector across all processes.

Warning
This method may only be called for locally replicated MultiVectors.
Precondition
isDistributed() == false

Definition at line 3554 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false > & Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::operator= ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type, false > &  source)

Assign the contents of source to this multivector (deep copy).

Precondition
The two multivectors must have the same communicator.
The input multivector's Map must be compatible with this multivector's Map. That is,
this->getMap ()->isCompatible (source.getMap ());
The two multivectors must have the same number of columns.
Note
This method must always be called as a collective operation on all processes over which the multivector is distributed. This is because the method reserves the right to check for compatibility of the two Maps, at least in debug mode.

Definition at line 2686 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subCopy ( const Teuchos::Range1D &  colRng) const

Return a MultiVector with copies of selected columns.

These methods are used to get the data underlying the MultiVector. They return data in one of three forms:

  • a MultiVector with a subset of the columns of the target MultiVector
  • a raw C pointer or array of raw C pointers
  • one of the Teuchos memory management classes 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 has constant stride.

Definition at line 2747 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subCopy ( const Teuchos::ArrayView< const size_t > &  cols) const

Return a MultiVector with copies of selected columns.

Definition at line 2716 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subView ( const Teuchos::Range1D &  colRng) const

Return a const MultiVector with const views of selected columns.

Definition at line 2954 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subView ( const Teuchos::ArrayView< const size_t > &  cols) const

Return a const MultiVector with const views of selected columns.

Definition at line 2906 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subViewNonConst ( const Teuchos::Range1D &  colRng)

Return a MultiVector with views of selected columns.

Definition at line 3075 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,node_type> > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::subViewNonConst ( const Teuchos::ArrayView< const size_t > &  cols)

Return a MultiVector with views of selected columns.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::offsetView ( const Teuchos::RCP< const map_type > &  subMap,
const size_t  offset 
) const

Return a const view of a subset of rows.

Return a const (nonmodifiable) view of this MultiVector consisting of a subset of the rows, as specified by an offset and a subset Map of this MultiVector's current row Map. If you want X1 or X2 to be nonconst (modifiable) views, use offsetViewNonConst() with the same arguments. "View" means "alias": if the original (this) MultiVector's data change, the view will see the changed data.

Parameters
subMap[in] The row Map for the new MultiVector. This must be a subset Map of this MultiVector's row Map.
offset[in] The local row offset at which to start the view.

Suppose that you have a MultiVector X, and you want to view X, on all processes in X's (MPI) communicator, as split into two row blocks X1 and X2. One could express this in Matlab notation as X = [X1; X2], except that here, X1 and X2 are views into X, rather than copies of X's data. This method assumes that the local indices of X1 and X2 are each contiguous, and that the local indices of X2 follow those of X1. If that is not the case, you cannot use views to divide X into blocks like this; you must instead use the Import or Export functionality, which copies the relevant rows of X.

Here is how you would construct the views X1 and X2.

using Teuchos::RCP;
MV X (...); // the input MultiVector
// ... fill X with data ...
// Map that on each process in X's communicator,
// contains the global indices of the rows of X1.
RCP<const map_type> map1 (new map_type (...));
// Map that on each process in X's communicator,
// contains the global indices of the rows of X2.
RCP<const map_type> map2 (new map_type (...));
// Create the first view X1. The second argument, the offset,
// is the index of the local row at which to start the view.
// X1 is the topmost block of X, so the offset is zero.
RCP<const MV> X1 = X.offsetView (map1, 0);
// Create the second view X2. X2 is directly below X1 in X,
// so the offset is the local number of rows in X1. This is
// the same as the local number of entries in map1.
RCP<const MV> X2 = X.offsetView (map2, X1->getLocalLength ());

It is legal, in the above example, for X1 or X2 to have zero local rows on any or all process(es). In that case, the corresponding Map must have zero local entries on that / those process(es). In particular, if X2 has zero local rows on a process, then the corresponding offset on that process would be the number of local rows in X (and therefore in X1) on that process. This is the only case in which the sum of the local number of entries in subMap (in this case, zero) and the offset may equal the number of local entries in *this.

Definition at line 2781 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::offsetViewNonConst ( const Teuchos::RCP< const map_type > &  subMap,
const size_t  offset 
)

Return a nonconst view of a subset of rows.

Return a nonconst (modifiable) view of this MultiVector consisting of a subset of the rows, as specified by an offset and a subset Map of this MultiVector's current row Map. If you want X1 or X2 to be const (nonmodifiable) views, use offsetView() with the same arguments. "View" means "alias": if the original (this) MultiVector's data change, the view will see the changed data, and if the view's data change, the original MultiVector will see the changed data.

Parameters
subMap[in] The row Map for the new MultiVector. This must be a subset Map of this MultiVector's row Map.
offset[in] The local row offset at which to start the view.

See the documentation of offsetView() for a code example and an explanation of edge cases.

Definition at line 2894 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getVector ( const size_t  j) const

Return a Vector which is a const view of column j.

Definition at line 3086 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getVectorNonConst ( const size_t  j)

Return a Vector which is a nonconst view of column j.

Definition at line 3113 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< const Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getData ( size_t  j) const

Const view of the local values in a particular vector of this multivector.

Definition at line 2596 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getDataNonConst ( size_t  j)

View of the local values in a particular vector of this multivector.

Definition at line 2639 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get1dCopy ( const Teuchos::ArrayView< Scalar > &  A,
const size_t  LDA 
) const

Fill the given array with a copy of this multivector's local values.

Parameters
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.

Definition at line 3124 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get2dCopy ( const Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > &  ArrayOfPtrs) const

Fill the given array with a copy of this multivector's local values.

Parameters
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.

Definition at line 3202 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< const Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get1dView ( ) const

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.

Definition at line 3240 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get2dView ( ) const

Return const persisting pointers to values.

Definition at line 3322 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get1dViewNonConst ( )

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.

Definition at line 3269 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::get2dViewNonConst ( )

Return non-const persisting pointers to values.

Definition at line 3298 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
KokkosClassic::MultiVector< Scalar, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getLocalMV ( ) const

A view of the underlying KokkosClassic::MultiVector object.

This method is for expert users only. It may change or be removed at any time.

Note
To Tpetra developers: This object's value type is Scalar and not impl_scalar_type, precisely because it supports a backwards compatibility use case.

Definition at line 3786 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dual_view_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getDualView ( ) const

Get the Kokkos::DualView which implements local storage.

Instead of getting the Kokkos::DualView, we highly recommend calling the templated view() method, that returns a Kokkos::View of the MultiVector's data in a given memory space. Since that MultiVector itself implements DualView semantics, it's much better to use MultiVector's interface to do "DualView things," like calling modify() and sync().

Warning
This method is ONLY for expert developers. Its interface may change or it may disappear at any time.

Definition at line 3825 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<class TargetDeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::sync ( )
inline

Update data on device or host only if data in the other space has been marked as modified.

If TargetDeviceType is the same as this MultiVector's device type, then copy data from host to device. Otherwise, copy data from device to host. In either case, only copy if the source of the copy has been modified.

This is a one-way synchronization only. If the target of the copy has been modified, this operation will discard those modifications. It will also reset both device and host modified flags.

Note
This method doesn't know on its own whether you modified the data in either memory space. You must manually mark the MultiVector as modified in the space in which you modified it, by calling the modify() method with the appropriate template parameter.

Definition at line 1139 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<class TargetDeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::modify ( )
inline

Mark data as modified on the given device TargetDeviceType.

If TargetDeviceType is the same as this MultiVector's device type, then mark the device's data as modified. Otherwise, mark the host's data as modified.

Definition at line 1149 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<class TargetDeviceType >
Kokkos::Impl::if_c< Kokkos::Impl::is_same< typename execution_space::memory_space, typename TargetDeviceType::memory_space>::value, typename dual_view_type::t_dev, typename dual_view_type::t_host>::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getLocalView ( ) const
inline

Return a view of the local data on a specific device.

Template Parameters
TargetDeviceTypeThe Kokkos Device type whose data to return.

Please don't be afraid of the if_c expression in the return value's type. That just tells the method what the return type should be: dual_view_type::t_dev if the TargetDeviceType template parameter matches this Tpetra object's device type, else dual_view_type::t_host.

For example, suppose you create a Tpetra::MultiVector for the Kokkos::Cuda device, like this:

typedef Kokkos::Compat::KokkosDeviceWrapperNode<Kokkos::Cuda> > node_type;
RCP<const map_type> map = ...;
mv_type DV (map, 3);

If you want to get the CUDA device Kokkos::View, do this:

typedef typename mv_type::dual_view_type dual_view_type;
typedef typename dual_view_type::t_dev device_view_type;
device_view_type cudaView = DV.getLocalView<Kokkos::Cuda> ();

and if you want to get the host mirror of that View, do this:

typedef typename dual_view_type::host_mirror_space host_execution_space;
typedef typename dual_view_type::t_host host_view_type;
host_view_type hostView = DV.getLocalView<host_execution_space> ();

Definition at line 1191 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Teuchos::ArrayView< dot_type > &  dots 
) const

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 $y^* x = \bar{y}^T x = \sum_i \bar{y}_i \cdot x_i$.

Precondition
*this and A have the same number of columns (vectors).
dots has at least as many entries as the number of columns in A.
Postcondition
dots[j] == (this->getVector[j])->dot (* (A.getVector[j]))

Definition at line 1302 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<dot_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Teuchos::ArrayView< T > &  dots 
) const
inline

Compute the dot product of each corresponding pair of vectors (columns) in A and B.

Template Parameters
TThe 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 1229 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<dot_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
std::vector< T > &  dots 
) const
inline

Like the above dot() overload, but for std::vector output.

Definition at line 1244 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Kokkos::View< dot_type *, execution_space > &  dots 
) const

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 $y^* x = \bar{y}^T x = \sum_i \bar{y}_i \cdot x_i$.

Parameters
A[in] MultiVector with which to dot *this.
dots[out] Device View with getNumVectors() entries.
Precondition
this->getNumVectors () == A.getNumVectors ()
dots.dimension_0 () == A.getNumVectors ()
Postcondition
dots(j) == (this->getVector[j])->dot (* (A.getVector[j]))

Definition at line 1203 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<dot_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Kokkos::View< T *, execution_space > &  dots 
) const
inline

Compute the dot product of each corresponding pair of vectors (columns) in A and B, storing the result in a device view.

Template Parameters
TThe 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 1291 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::abs ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A)

Put element-wise absolute values of input Multi-vector in target: A = abs(this)

Definition at line 2387 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::reciprocal ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A)

Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).

Definition at line 2331 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::scale ( const Scalar &  alpha)

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.

Definition at line 2082 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::scale ( const Teuchos::ArrayView< const Scalar > &  alpha)

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.

Definition at line 2143 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::scale ( const Kokkos::View< const impl_scalar_type *, execution_space > &  alpha)

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.

Definition at line 2169 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::scale ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A 
)

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.

Definition at line 2247 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::update ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Scalar &  beta 
)

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.

Definition at line 2433 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::update ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const Scalar &  beta,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  B,
const Scalar &  gamma 
)

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.

Definition at line 2518 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm1 ( const Kokkos::View< mag_type *, execution_space > &  norms) const

Compute the one-norm of each vector (column), storing the result in a device view.

The one-norm of a vector is the sum of squares of the magnitudes of the vector's entries. On exit, norms(j) is the one-norm of column j of this MultiVector.

Parameters
norms[out] Device View with getNumVectors() entries.
Precondition
norms.dimension_0 () == this->getNumVectors ()
Postcondition
norms(j) == (this->getVector[j])->norm1 (* (A.getVector[j]))

Definition at line 1457 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm1 ( const Kokkos::View< T *, execution_space > &  norms) const
inline

Compute the one-norm of each vector (column), storing the result in a device view.

Template Parameters
TThe output type of the dot products.

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 1406 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm1 ( const Teuchos::ArrayView< mag_type > &  norms) const

Compute the one-norm of each vector (column).

The one-norm of a vector is the sum of squares of the magnitudes of the vector's entries. On exit, norms[j] is the one-norm of column j of this MultiVector.

Definition at line 1439 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type,T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm1 ( const Teuchos::ArrayView< T > &  norms) const
inline

Compute the one-norm of each vector (column).

Template Parameters
TThe output type of the norms.

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 1442 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm2 ( const Kokkos::View< mag_type *, execution_space > &  norms) const

Compute the two-norm of each vector (column), storing the result in a device view.

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.

Parameters
norms[out] Device View with getNumVectors() entries.
Precondition
norms.dimension_0 () == this->getNumVectors ()
Postcondition
norms(j) == (this->getVector[j])->dot (* (A.getVector[j]))

Definition at line 1338 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm2 ( const Kokkos::View< T *, execution_space > &  norms) const
inline

Compute the two-norm of each vector (column), storing the result in a device view.

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 1481 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm2 ( const Teuchos::ArrayView< mag_type > &  norms) const

Compute the two-norm of each vector (column).

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.

Definition at line 1320 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type,T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::norm2 ( const Teuchos::ArrayView< T > &  norms) const
inline

Compute the two-norm of each vector (column).

Template Parameters
TThe output type of the norms.

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 1517 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normInf ( const Kokkos::View< mag_type *, execution_space > &  norms) const

Compute the infinity-norm of each vector (column), storing the result in a device 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.

Definition at line 1484 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type, T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normInf ( const Kokkos::View< T *, execution_space > &  norms) const
inline

Compute the two-norm of each vector (column), storing the result in a device view.

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 1550 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normInf ( const Teuchos::ArrayView< mag_type > &  norms) const

Compute the infinity-norm of each vector (column).

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.

Definition at line 1466 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type,T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normInf ( const Teuchos::ArrayView< T > &  norms) const
inline

Compute the infinity-norm of each vector (column).

Template Parameters
TThe output type of the norms.

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 1585 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normWeighted ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  weights,
const Teuchos::ArrayView< mag_type > &  norms 
) const

Compute Weighted 2-norm (RMS Norm) of each vector in multi-vector. The outcome of this routine is undefined for non-floating point scalar types (e.g., int).

Definition at line 1347 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<typename T >
Kokkos::Impl::enable_if< !(Kokkos::Impl::is_same<mag_type,T>::value), void >::type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normWeighted ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  weights,
const Teuchos::ArrayView< T > &  norms 
) const
inline

Compute the weighted 2-norm (RMS Norm) of each column.

The outcome of this routine is undefined for non-floating point scalar types (e.g., int).

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 1619 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::meanValue ( const Teuchos::ArrayView< impl_scalar_type > &  means) const

Compute mean (average) value of each column.

The outcome of this routine is undefined for non-floating point scalar types (e.g., int).

Definition at line 1734 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::multiply ( Teuchos::ETransp  transA,
Teuchos::ETransp  transB,
const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  B,
const Scalar &  beta 
)

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.

Definition at line 3346 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::elementWiseMultiply ( Scalar  scalarAB,
const Vector< Scalar, LocalOrdinal, GlobalOrdinal, node_type, false > &  A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  B,
Scalar  scalarThis 
)

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:

C(i,j) = scalarThis * C(i,j) + scalarAB * B(i,j) * A(i,1);

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.

Definition at line 3497 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getNumVectors ( ) const

Number of columns in the multivector.

Definition at line 1060 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getLocalLength ( ) const

Local number of rows on the calling process.

Definition at line 664 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
global_size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getGlobalLength ( ) const

Global number of rows in the multivector.

Definition at line 676 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getStride ( ) const

Stride between columns in the multivector.

This is only meaningful if isConstantStride() returns true.

Warning
This may be different on different processes.

Definition at line 688 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::isConstantStride ( ) const

Whether this multivector has constant stride between columns.

Warning
This may be different on different processes.

Definition at line 657 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
std::string Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::description ( ) const
virtual

A simple one-line description of this object.

Definition at line 3834 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
virtual

Print the object with the given verbosity level to a FancyOStream.

Parameters
out[out] Output stream to which to print. For verbosity levels VERB_LOW and lower, only the process with rank 0 ("Proc 0") in the MultiVector's communicator prints. For verbosity levels strictly higher than VERB_LOW, all processes in the communicator need to be able to print to the output stream.
verbLevel[in] Verbosity level. The default verbosity (verbLevel=VERB_DEFAULT) is VERB_LOW.

The amount and content of what this method prints depends on the verbosity level. In the list below, each higher level includes all the content of the previous levels, as well as its own content.

  • VERB_LOW: Only Proc 0 prints; it prints the same thing as description().
  • VERB_MEDIUM: Each process prints its local length (the number of rows that it owns).
  • VERB_HIGH: Each process prints whether the multivector has constant stride (see isConstantStride()), and if so, what that stride is. (Stride may differ on different processes.)
  • VERB_EXTREME: Each process prints the values in its local part of the multivector. This will print out as many rows of data as the global number of rows in the multivector, so beware.

Definition at line 3855 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, node_type > > &  newMap)
virtual

Remove processes owning zero rows from the Map and their communicator.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in Tpetra_DistObject_decl.hpp.
We make NO promises of backwards compatibility. This method may change or disappear at any time.
Parameters
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.

Definition at line 4009 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::setCopyOrView ( const Teuchos::DataAccess  copyOrView)
inline

Set whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.

Warning
The Kokkos refactor version of MultiVector only implements view semantics. If you attempt to call this method with copyOrView == Teuchos::Copy, it will throw std::invalid_argument.
This method is only for expert use. It may change or disappear at any time.

Definition at line 1783 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::DataAccess Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getCopyOrView ( ) const
inline

Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.

Note
The Kokkos refactor version of MultiVector only implements view semantics. This is not currently true for the "classic" version of MultiVector, though that will change in the near future.
Warning
This method is only for expert use. It may change or disappear at any time.

Definition at line 1802 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::assign ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, node_type > &  src)

Copy the contents of src into *this (deep copy).

Parameters
src[in] Source MultiVector (input of the deep copy).
Precondition
! src.getMap ().is_null () && ! this->getMap ().is_null ()
src.getMap ()->isCompatible (* (this->getMap ())
Postcondition
Any outstanding views of src or *this remain valid.
Note
To implementers: The postcondition implies that the implementation must not reallocate any memory of *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.

Definition at line 4019 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::normImpl ( const Kokkos::View< mag_type *, execution_space > &  norms,
const EWhichNorm  whichNorm 
) const
protected

Compute the norm of each vector (column), storing the result in a device View.

This method consolidates all common code between the infinity-norm, 1-norm, and 2-norm calculations. On exit, norms(j) is the norm (of the selected type) of column j of this MultiVector.

Definition at line 1664 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
template<class T >
Teuchos::ArrayRCP< T > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getSubArrayRCP ( Teuchos::ArrayRCP< T >  arr,
size_t  j 
) const
protected

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.

Definition at line 3765 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getOrigNumLocalRows ( ) const
protected

"Original" number of rows in the (local) data.

Definition at line 2763 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::getOrigNumLocalCols ( ) const
protected

"Original" number of columns in the (local) data.

Definition at line 2772 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::checkSizes ( const SrcDistObject sourceObj)
protectedvirtual

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.

Definition at line 706 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::constantNumberOfPackets ( ) const
protectedvirtual

Number of packets to send per LID.

Definition at line 729 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
virtual bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::useNewInterface ( )
inlineprotectedvirtual

Whether this class implements the old or new interface of DistObject.

Definition at line 1943 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

Member Data Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
dual_view_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::view_
mutableprotected

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 1839 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
dual_view_type Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::origView_
mutableprotected

The "original view" of the MultiVector's data.

Methods like offsetView() return a view of a contiguous subset of rows. At some point, we might like to get all of the rows back, by taking another view of a superset of rows. For example, we might like to get a column Map view of a (domain Map view of a (column Map MultiVector)). Tpetra's implementation of Gauss-Seidel and SOR in CrsMatrix relies on this functionality. However, Kokkos (rightfully) forbids us from taking a superset of rows of the current view.

We deal with this at the Tpetra level by keeping around the original view of all the rows (and columns), which is origView_. Methods like offsetView() then use origView_, not view_, to make the subview for the returned MultiVector. Furthermore, offsetView() can do error checking by getting the original number of rows from origView_.

This may pose some problems for offsetView if it is given an offset other than zero, but that case is hardly exercised, so I am not going to worry about it for now.

Note that the "original" view isn't always original. It always has the original number of rows. However, some special cases of constructors that take a whichVectors argument, when whichVectors.size() is 1, may point origView_ to the column to view. Those constructors do this so that the resulting MultiVector has constant stride. This special case does not affect correctness of offsetView and related methods.

Definition at line 1870 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::Array<size_t> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType >, false >::whichVectors_
protected

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 1884 of file Tpetra_KokkosRefactor_MultiVector_decl.hpp.


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