Tpetra Matrix/Vector Services  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines
Public Member Functions | Protected Member Functions | Protected Attributes | Related Functions
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

A distributed dense vector. More...

#include <Tpetra_Vector_decl.hpp>

Inheritance diagram for Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
Inheritance graph
[legend]

List of all members.

Public Types

Typedefs to facilitate template metaprogramming
typedef Scalar scalar_type
 The type of each entry in the Vector.
typedef LocalOrdinal local_ordinal_type
 The type of local indices.
typedef GlobalOrdinal global_ordinal_type
 The type of global indices.
typedef Node node_type
 The Kokkos Node type.
typedef base_type::dot_type dot_type
 Type of an inner ("dot") product result.
typedef base_type::mag_type mag_type
 Type of a norm result.
typedef Map
< local_ordinal_type,
global_ordinal_type, node_type
map_type
 The type of the Map specialization used by this class.
Typedefs to facilitate template metaprogramming.
typedef DistObject< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > 
DO
Typedefs
typedef Scalar packet_type
 The type of each datum being sent or received in an Import or Export.

Public Member Functions

virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
 Remove processes owning zero rows from the Map and their communicator.
void setCopyOrView (const Teuchos::DataAccess copyOrView)
 Set whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.
Teuchos::DataAccess getCopyOrView () const
 Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics.
void assign (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &src)
 Copy the contents of src into *this (deep copy).
Constructors and destructor
 Vector (const Teuchos::RCP< const map_type > &map, const bool zeroOut=true)
 Basic constructor.
 Vector (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
 Copy constructor.
 Vector (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::DataAccess copyOrView)
 Copy constructor, with option to do shallow copy and mark the result as having "view semantics.".
 Vector (const Teuchos::RCP< const map_type > &map, const Teuchos::ArrayView< const Scalar > &A)
 Set vector values from an existing array (copy)
virtual ~Vector ()
 Destructor.
Clone method
template<class Node2 >
Teuchos::RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node2 > > 
clone (const Teuchos::RCP< Node2 > &node2)
 Return a deep copy of *this with a different Node type.
Post-construction modification routines
void replaceGlobalValue (GlobalOrdinal globalRow, const Scalar &value)
 Replace current value at the specified location with specified value.
void sumIntoGlobalValue (GlobalOrdinal globalRow, const Scalar &value)
 Adds specified value to existing value at the specified location.
void replaceLocalValue (LocalOrdinal myRow, const Scalar &value)
 Replace current value at the specified location with specified values.
void sumIntoLocalValue (LocalOrdinal myRow, const Scalar &value)
 Adds specified value to existing value at the specified location.
Extraction methods
void get1dCopy (Teuchos::ArrayView< Scalar > A) const
 Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management classes).
Teuchos::ArrayRCP< Scalar > getDataNonConst ()
 View of the local values of this vector.
Teuchos::ArrayRCP< const Scalar > getData () const
 Const view of the local values of this vector.
Mathematical methods
Scalar dot (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &a) const
 Computes dot product of this Vector against input Vector x.
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm1 () const
 Return 1-norm of this Vector.
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm2 () const
 Compute 2-norm of this Vector.
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
normInf () const
 Compute Inf-norm of this Vector.
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
normWeighted (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights) const
 Compute Weighted 2-norm (RMS Norm) of this Vector.
Scalar meanValue () const
 Compute mean (average) value of this Vector.
Methods to get a row subset view of a Vector (as a Vector).

These methods hide those of the same names in MultiVector. This is necessary because users who use Vectors would very much like a view of a Vector to be a Vector, rather than a MultiVector. Not only do these methods in MultiVector return MultiVector pointers, those pointers can't even be dynamic_cast to Vector, because the returned objects are not Vectors.

There are two options to make this work:

  1. Make these methods virtual, and override them in Vector.
  2. Hide these methods in MultiVector (that return MultiVector) with methods in Vector that return Vector.

The first option only solves the dynamic_cast problem mentioned above. The methods still would have to return MultiVector pointers. Thus, we prefer the second option. Even though the methods in Vector must therefore hide those in MultiVector, Vector is a MultiVector, so the returned objects from the Vector methods are still compatible with MultiVector operations.

Teuchos::RCP< const Vector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
offsetView (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &subMap, size_t offset) const
 Return a const view of a subset of rows.
Teuchos::RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
offsetViewNonConst (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &subMap, size_t offset)
 Return a nonconst view of a subset of rows.
Implementation of the Teuchos::Describable interface
virtual std::string description () const
 A simple one-line description of this object.
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print the object with some verbosity level to a FancyOStream.
Constructors and destructor
template<class Node2 >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node2 > > 
clone (const Teuchos::RCP< Node2 > &node2) const
 Create a cloned MultiVector for a different node type.
Post-construction modification routines
void replaceGlobalValue (GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
 Replace value, using global (row) index.
void sumIntoGlobalValue (GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
 Add value to existing value, using global (row) index.
void replaceLocalValue (LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
 Replace value, using local (row) index.
void sumIntoLocalValue (LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
 Add value to existing value, using local (row) index.
void putScalar (const Scalar &value)
 Set all values in the multivector with the given value.
void randomize ()
 Set all values in the multivector to pseudorandom numbers.
void replaceMap (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
 Replace the underlying Map in place.
void reduce ()
 Sum values of a locally replicated multivector across all processes.
Data copy and view methods

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

  • A MultiVector with a subset of the rows or columns of the input MultiVector
  • An array of data, wrapped in 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.

Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
subCopy (const Teuchos::Range1D &colRng) const
 Return a MultiVector with copies of selected columns.
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
subCopy (const Teuchos::ArrayView< const size_t > &cols) const
 Return a MultiVector with copies of selected columns.
Teuchos::RCP< const
MultiVector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
subView (const Teuchos::Range1D &colRng) const
 Return a MultiVector with const views of selected columns.
Teuchos::RCP< const
MultiVector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
subView (const Teuchos::ArrayView< const size_t > &cols) const
 Return a const MultiVector with const views of selected columns.
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
subViewNonConst (const Teuchos::Range1D &colRng)
 Return a MultiVector with views of selected columns.
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
subViewNonConst (const Teuchos::ArrayView< const size_t > &cols)
 Return a MultiVector with views of selected columns.
Teuchos::RCP< const Vector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
getVector (size_t j) const
 Return a Vector which is a const view of column j.
Teuchos::RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
getVectorNonConst (size_t j)
 Return a Vector which is a nonconst view of column j.
Teuchos::ArrayRCP< const Scalar > getData (size_t j) const
 Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > getDataNonConst (size_t j)
 View of the local values in a particular vector of this multivector.
void get1dCopy (Teuchos::ArrayView< Scalar > A, size_t LDA) const
 Fill the given array with a copy of this multivector's local values.
void get2dCopy (Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
 Fill the given array with a copy of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > get1dView () const
 Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP
< Teuchos::ArrayRCP< const
Scalar > > 
get2dView () const
 Return const persisting pointers to values.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst ()
 Nonconst persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP
< Teuchos::ArrayRCP< Scalar > > 
get2dViewNonConst ()
 Return non-const persisting pointers to values.
KokkosClassic::MultiVector
< Scalar, Node > 
getLocalMV () const
 A view of the underlying KokkosClassic::MultiVector object.
TEUCHOS_DEPRECATED
KokkosClassic::MultiVector
< Scalar, Node > & 
getLocalMVNonConst ()
 A nonconst reference to a view of the underlying KokkosClassic::MultiVector object.
Mathematical methods
void dot (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
 Compute the dot product of each corresponding pair of vectors (columns) in A and B.
void abs (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
 Put element-wise absolute values of input Multi-vector in target: A = abs(this)
void reciprocal (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
 Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
void scale (const Scalar &alpha)
 Scale in place: this = alpha*this.
void scale (Teuchos::ArrayView< const Scalar > alpha)
 Scale each column in place: this[j] = alpha[j]*this[j].
void scale (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
 Scale in place: this = alpha * A.
void update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
 Update: this = beta*this + alpha*A.
void update (const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma)
 Update: this = gamma*this + alpha*A + beta*B.
void norm1 (const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
 Compute the one-norm of each vector (column).
void norm2 (const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
 Compute the two-norm of each vector (column).
void normInf (const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
 Compute the infinity-norm of each vector (column).
void normWeighted (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &weights, const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void meanValue (const Teuchos::ArrayView< Scalar > &means) const
 Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined for non-floating point scalar types (e.g., int).
void multiply (Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta)
 Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
void elementWiseMultiply (Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
 Multiply a Vector A elementwise by a MultiVector B.
Attribute access functions
size_t getNumVectors () const
 Number of columns in the multivector.
size_t getLocalLength () const
 Local number of rows on the calling process.
global_size_t getGlobalLength () const
 Global number of rows in the multivector.
size_t getStride () const
 Stride between columns in the multivector.
bool isConstantStride () const
 Whether this multivector has constant stride between columns.
Public methods for redistributing data
void doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 Import data into this object using an Import object ("forward mode").
void doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Import data into this object using an Export object ("reverse mode").
void doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
 Export data into this object using an Export object ("forward mode").
void doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
 Export data into this object using an Import object ("reverse mode").
Attribute accessor methods
bool isDistributed () const
 Whether this is a globally distributed object.
virtual Teuchos::RCP< const
Map< LocalOrdinal,
GlobalOrdinal, Node > > 
getMap () const
 The Map describing the parallel distribution of this object.
I/O methods
void print (std::ostream &os) const
 Print this object to the given output stream.

Protected Member Functions

 Vector (const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const ArrayRCP< Scalar > &data)
 Advanced constructor accepting parallel buffer view, used by MultiVector to break off Vector objects.
 Vector (const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const ArrayRCP< Scalar > &data, EPrivateComputeViewConstructor)
 Advanced constructor accepting parallel buffer view, used by MultiVector to break off Vector objects.
virtual void doTransfer (const SrcDistObject &src, CombineMode CM, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs, const Teuchos::ArrayView< const LocalOrdinal > &remoteLIDs, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Distributor &distor, ReverseOption revOp)
 Redistribute data across memory images.
Methods implemented by subclasses and used by doTransfer().

The doTransfer() method uses the subclass' implementations of these methods to implement data transfer. Subclasses of DistObject must implement these methods. This is an instance of the Template Method Pattern. ("Template" here doesn't mean "C++ template"; it means "pattern with holes that are filled in by the subclass' method implementations.")

virtual void copyAndPermute (const SrcDistObject &source, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs)=0
 Perform copies and permutations that are local to this process.
virtual void packAndPrepare (const SrcDistObject &source, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< Scalar > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor)=0
 Perform any packing or preparation required for communication.
virtual void unpackAndCombine (const Teuchos::ArrayView< const LocalOrdinal > &importLIDs, const Teuchos::ArrayView< const Scalar > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM)=0
 Perform any unpacking and combining after communication.

Protected Attributes

KMV lclMV_
 The KokkosClassic::MultiVector containing the compute buffer of data.
Array< size_t > whichVectors_
 Indices of columns this multivector is viewing.
bool hasViewSemantics_
 Whether this MultiVector has view semantics.
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
map_
 The Map over which this object is distributed.
Teuchos::Array< Scalar > imports_
 Buffer into which packed data are imported (received from other processes).
Teuchos::Array< size_t > numImportPacketsPerLID_
 Number of packets to receive for each receive operation.
Teuchos::Array< Scalar > exports_
 Buffer from which packed data are exported (sent to other processes).
Teuchos::Array< size_t > numExportPacketsPerLID_
 Number of packets to send for each send operation.

Friends

Methods for use only by experts
void removeEmptyProcessesInPlace (Teuchos::RCP< Tpetra::DistObject< PT, LO, GO, NT > > &input, const Teuchos::RCP< const Map< LO, GO, NT > > &newMap)
void removeEmptyProcessesInPlace (Teuchos::RCP< Tpetra::DistObject< PT, LO, GO, NT > > &input)

Related Functions

(Note that these are not member functions.)

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal,
Kokkos::Compat::KokkosDeviceWrapperNode
< DeviceType > > > 
createMultiVectorFromView (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > &map, const Teuchos::ArrayRCP< Scalar > &view, const size_t LDA, const size_t numVectors)
 Nonmember MultiVector constructor with view semantics using user-allocated data.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Kokkos::Compat::KokkosDeviceWrapperNode
< DeviceType > > > 
createVectorFromView (const RCP< const Map< LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > &map, const ArrayRCP< Scalar > &view)
 Non-member function to create a Vector from a specified Map.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
createMultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t numVectors)
 Nonmember MultiVector constructor: make a MultiVector from a given Map.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
createMultiVectorFromView (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayRCP< Scalar > &view, const size_t LDA, const size_t numVectors)
 Nonmember MultiVector constructor with view semantics using user-allocated data.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
createVector (const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
 Non-member function to create a Vector from a specified Map.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
createVectorFromView (const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const ArrayRCP< Scalar > &view)
 Non-member function to create a Vector with view semantics using user-allocated data.
template<class DS , class DL , class DG , class DN , class SS , class SL , class SG , class SN >
void deep_copy (MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
 Copy the contents of the MultiVector src into dst.
template<class ST , class LO , class GO , class NT >
MultiVector< ST, LO, GO, NT > createCopy (const MultiVector< ST, LO, GO, NT > &src)
 Return a deep copy of the MultiVector src.

View constructors, used only by nonmember constructors.

struct Details::CreateMultiVectorFromView
bool vectorIndexOutOfRange (size_t VectorIndex) const
template<class T >
ArrayRCP< T > getSubArrayRCP (ArrayRCP< T > arr, size_t j) const
 Persisting view of j-th column in the given ArrayRCP.

Implementation of Tpetra::DistObject

virtual bool checkSizes (const SrcDistObject &sourceObj)
 Whether data redistribution between sourceObj and this object is legal.
virtual size_t constantNumberOfPackets () const
 Number of packets to send per LID.
virtual void copyAndPermute (const SrcDistObject &sourceObj, size_t numSameIDs, const ArrayView< const LocalOrdinal > &permuteToLIDs, const ArrayView< const LocalOrdinal > &permuteFromLIDs)
virtual void packAndPrepare (const SrcDistObject &sourceObj, const ArrayView< const LocalOrdinal > &exportLIDs, Array< Scalar > &exports, const ArrayView< size_t > &numExportPacketsPerLID, size_t &constantNumPackets, Distributor &distor)
virtual void unpackAndCombine (const ArrayView< const LocalOrdinal > &importLIDs, const ArrayView< const Scalar > &imports, const ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM)
void createViews () const
 Hook for creating a const view.
void createViewsNonConst (KokkosClassic::ReadWriteOption rwo)
 Hook for creating a nonconst view.
void releaseViews () const
 Hook for releasing views.
ArrayRCP< Scalar > ncview_
 Nonconst host view created in createViewsNonConst().
ArrayRCP< const Scalar > cview_
 Const host view created in createViews().

Detailed Description

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
class Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >

A distributed dense vector.

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

This class inherits from MultiVector, and has the same template parameters. A Vector is a special case of a MultiVector that has only one vector (column). It may be used wherever a MultiVector may be used. Please see the documentation of MultiVector for more details.

Examples:
lesson02_init_map_vec.cpp, lesson02_read_modify_vec.cpp, and lesson03_power_method.cpp.

Definition at line 72 of file Tpetra_Vector_decl.hpp.


Member Typedef Documentation

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef Scalar Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type

The type of each entry in the Vector.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 85 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef LocalOrdinal Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type

The type of local indices.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 87 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef GlobalOrdinal Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type

The type of global indices.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 89 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef Node Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type

The Kokkos Node type.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 91 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef base_type::dot_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot_type

Type of an inner ("dot") product result.

This is usually the same as scalar_type, but may differ if scalar_type is e.g., an uncertainty quantification type from the Stokhos package.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 98 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef base_type::mag_type Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type

Type of a norm result.

This is usually the same as the type of the magnitude (absolute value) of scalar_type, but may differ if scalar_type is e.g., an uncertainty quantification type from the Stokhos package.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 106 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
typedef Map<local_ordinal_type, global_ordinal_type, node_type> Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::map_type

The type of the Map specialization used by this class.

Definition at line 109 of file Tpetra_Vector_decl.hpp.

typedef Scalar Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::packet_type [inherited]

The type of each datum being sent or received in an Import or Export.

Note that this type does not always correspond to the Scalar template parameter of subclasses.

Definition at line 183 of file Tpetra_DistObject_decl.hpp.


Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const Teuchos::RCP< const map_type > &  map,
const bool  zeroOut = true 
) [explicit]

Basic constructor.

Parameters:
map[in] The Vector's Map. The Map describes the distribution of rows over process(es) in the Map's communicator.
zeroOut[in] If true (the default), require that all the Vector's entries be zero on return. If false, the Vector's entries have undefined values on return, and must be set explicitly.

Definition at line 57 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source)

Copy constructor.

Whether this does a deep copy or a shallow copy depends on whether source has "view semantics." See discussion in the documentation of the two-argument copy constructor below.

Definition at line 63 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Teuchos::DataAccess  copyOrView 
)

Copy constructor, with option to do shallow copy and mark the result as having "view semantics.".

If copyOrView is Teuchos::View, this constructor marks the result as having "view semantics." This means that copy construction or assignment (operator=) with the resulting object will always do a shallow copy, and will transmit view semantics to the result of the shallow copy. If copyOrView is Teuchos::Copy, this constructor always does a deep copy and marks the result as not having view semantics, whether or not source has view semantics.

View semantics are a "forwards compatibility" measure for porting to the Kokkos refactor version of Tpetra. The latter only ever has view semantics. The "classic" version of Tpetra does not currently have view semantics by default, but this will change.

Definition at line 69 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const Teuchos::RCP< const map_type > &  map,
const Teuchos::ArrayView< const Scalar > &  A 
)

Set vector values from an existing array (copy)

Definition at line 92 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~Vector ( ) [virtual]

Destructor.

Definition at line 99 of file Tpetra_Vector_def.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const ArrayRCP< Scalar > &  data 
) [protected]

Advanced constructor accepting parallel buffer view, used by MultiVector to break off Vector objects.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Vector ( const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const ArrayRCP< Scalar > &  data,
EPrivateComputeViewConstructor   
) [protected]

Advanced constructor accepting parallel buffer view, used by MultiVector to break off Vector objects.


Member Function Documentation

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
template<class Node2 >
Teuchos::RCP<Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::clone ( const Teuchos::RCP< Node2 > &  node2)

Return a deep copy of *this with a different Node type.

Template Parameters:
Node2The returned Vector's Node type.
Parameters:
node2[in] The returned Vector's Kokkos Node instance.
template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValue ( GlobalOrdinal  globalRow,
const Scalar &  value 
)

Replace current value at the specified location with specified value.

Precondition:
globalRow must be a valid global element on this node, according to the row map.

Definition at line 103 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValue ( GlobalOrdinal  globalRow,
const Scalar &  value 
)

Adds specified value to existing value at the specified location.

Precondition:
globalRow must be a valid global element on this node, according to the row map.

Definition at line 108 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValue ( LocalOrdinal  myRow,
const Scalar &  value 
)

Replace current value at the specified location with specified values.

Precondition:
localRow must be a valid local element on this node, according to the row map.

Definition at line 113 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValue ( LocalOrdinal  myRow,
const Scalar &  value 
)

Adds specified value to existing value at the specified location.

Precondition:
localRow must be a valid local element on this node, according to the row map.

Definition at line 118 of file Tpetra_Vector_def.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get1dCopy ( Teuchos::ArrayView< Scalar >  A) const

Return multi-vector values in user-provided two-dimensional array (using Teuchos memory management classes).

Definition at line 123 of file Tpetra_Vector_def.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
Teuchos::ArrayRCP<Scalar> Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDataNonConst ( ) [inline]

View of the local values of this vector.

Definition at line 209 of file Tpetra_Vector_decl.hpp.

template<class Scalar = MultiVector<>::scalar_type, class LocalOrdinal = typename MultiVector<Scalar>::local_ordinal_type, class GlobalOrdinal = typename MultiVector<Scalar, LocalOrdinal>::global_ordinal_type, class Node = typename MultiVector<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
Teuchos::ArrayRCP<const Scalar> Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getData ( ) const [inline]

Const view of the local values of this vector.

Definition at line 215 of file Tpetra_Vector_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Scalar Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  a) const

Computes dot product of this Vector against input Vector x.

Definition at line 129 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ScalarTraits< Scalar >::magnitudeType Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm1 ( ) const

Return 1-norm of this Vector.

Definition at line 161 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ScalarTraits< Scalar >::magnitudeType Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm2 ( ) const

Compute 2-norm of this Vector.

Definition at line 170 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ScalarTraits< Scalar >::magnitudeType Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normInf ( ) const

Compute Inf-norm of this Vector.

Definition at line 179 of file Tpetra_Vector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::ScalarTraits< Scalar >::magnitudeType Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normWeighted ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  weights) const

Compute Weighted 2-norm (RMS Norm) of this Vector.

Definition at line 189 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Scalar Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::meanValue ( ) const

Compute mean (average) value of this Vector.

Definition at line 151 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::offsetView ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  subMap,
size_t  offset 
) const

Return a const view of a subset of rows.

Return a const (nonmodifiable) view of this Vector consisting of a subset of the rows, as specified by an offset and a subset Map of this Vector'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 Vector. 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 Vector 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;
 typedef Tpetra::Map<LO, GO, Node> map_type;
 typedef Tpetra::Vector<Scalar, LO, GO, Node> V;

 V X (...); // the input Vector
 // ... 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 V> 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 V> 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.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 200 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::offsetViewNonConst ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  subMap,
size_t  offset 
)

Return a nonconst view of a subset of rows.

Return a nonconst (modifiable) view of this Vector consisting of a subset of the rows, as specified by an offset and a subset Map of this Vector'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) Vector's data change, the view will see the changed data, and if the view's data change, the original Vector will see the changed data.

Parameters:
subMap[in] The row Map for the new Vector. This must be a subset Map of this Vector'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.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 229 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::string Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const [virtual]

A simple one-line description of this object.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 257 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const [virtual]

Print the object with some verbosity level to a FancyOStream.

Reimplemented from Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 279 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
template<class Node2 >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node2 > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::clone ( const Teuchos::RCP< Node2 > &  node2) const [inherited]

Create a cloned MultiVector for a different node type.

Definition at line 1524 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValue ( GlobalOrdinal  globalRow,
size_t  vectorIndex,
const Scalar &  value 
) [inherited]

Replace value, using global (row) index.

Replace the current value at row globalRow (a global index) and column vectorIndex with the given value. The column index is zero based.

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

Definition at line 3021 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValue ( GlobalOrdinal  globalRow,
size_t  vectorIndex,
const Scalar &  value 
) [inherited]

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 vectorIndex. The column index is zero based.

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

Definition at line 3045 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValue ( LocalOrdinal  myRow,
size_t  vectorIndex,
const Scalar &  value 
) [inherited]

Replace value, using local (row) index.

Replace the current value at row myRow (a local index) and column vectorIndex with the given value. The column index is zero based.

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

Definition at line 2963 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValue ( LocalOrdinal  myRow,
size_t  vectorIndex,
const Scalar &  value 
) [inherited]

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

Add the given value to the existing value at row myRow (a local index) and column vectorIndex. The column index is zero based.

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

Definition at line 2992 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::putScalar ( const Scalar &  value) [inherited]

Set all values in the multivector with the given value.

Definition at line 1592 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::randomize ( ) [inherited]

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.
Examples:
lesson03_power_method.cpp.

Definition at line 1571 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceMap ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map) [inherited]

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 1613 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::reduce ( ) [inherited]

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 2896 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subCopy ( const Teuchos::Range1D &  colRng) const [inherited]

Return a MultiVector with copies of selected columns.

Parameters:
colRng[in] Inclusive, contiguous range of columns. [colRng.lbound(), colRng.ubound()] defines the range.

Definition at line 2175 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subCopy ( const Teuchos::ArrayView< const size_t > &  cols) const [inherited]

Return a MultiVector with copies of selected columns.

Definition at line 2125 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subView ( const Teuchos::Range1D &  colRng) const [inherited]

Return a MultiVector with const views of selected columns.

Parameters:
colRng[in] Inclusive, contiguous range of columns. [colRng.lbound(), colRng.ubound()] defines the range.

Definition at line 2366 of file Tpetra_MultiVector_def.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
Teuchos::RCP<const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subView ( const Teuchos::ArrayView< const size_t > &  cols) const [inherited]

Return a const MultiVector with const views of selected columns.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subViewNonConst ( const Teuchos::Range1D &  colRng) [inherited]

Return a MultiVector with views of selected columns.

Parameters:
colRng[in] Inclusive, contiguous range of columns. [colRng.lbound(), colRng.ubound()] defines the range.

Definition at line 2455 of file Tpetra_MultiVector_def.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::subViewNonConst ( const Teuchos::ArrayView< const size_t > &  cols) [inherited]

Return a MultiVector with views of selected columns.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getVector ( size_t  j) const [inherited]

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

Definition at line 2487 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getVectorNonConst ( size_t  j) [inherited]

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

Definition at line 2512 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< const Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getData ( size_t  j) const [inherited]

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

Definition at line 2027 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDataNonConst ( size_t  j) [inherited]

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

Definition at line 2037 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get1dCopy ( Teuchos::ArrayView< Scalar >  A,
size_t  LDA 
) const [inherited]

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 2533 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get2dCopy ( Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > >  ArrayOfPtrs) const [inherited]

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 2569 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< const Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get1dView ( ) const [inherited]

Const persisting (1-D) view of this multivector's local values.

This method assumes that the columns of the multivector are stored contiguously. If not, this method throws std::runtime_error.

Definition at line 2604 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get2dView ( ) const [inherited]

Return const persisting pointers to values.

Definition at line 2685 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< Scalar > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get1dViewNonConst ( ) [inherited]

Nonconst persisting (1-D) view of this multivector's local values.

This method assumes that the columns of the multivector are stored contiguously. If not, this method throws std::runtime_error.

Definition at line 2625 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::get2dViewNonConst ( ) [inherited]

Return non-const persisting pointers to values.

Definition at line 2647 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
KokkosClassic::MultiVector< Scalar, Node > Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMV ( ) const [inherited]

A view of the underlying KokkosClassic::MultiVector object.

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

Definition at line 3087 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
TEUCHOS_DEPRECATED KokkosClassic::MultiVector< Scalar, Node > & Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMVNonConst ( ) [inherited]

A nonconst reference to a view of the underlying KokkosClassic::MultiVector object.

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

Warning:
This method is DEPRECATED. It may disappear at any time. Please call getLocalMV() instead. There was never actually a need for a getLocalMVNonConst() method, as far as I can tell.

Definition at line 3094 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dot ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Teuchos::ArrayView< Scalar > &  dots 
) const [inherited]

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

The "dot product" is the standard Euclidean inner product. If the type of entries of the vectors (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 1244 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::abs ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A) [inherited]

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

Definition at line 1896 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::reciprocal ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A) [inherited]

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

Definition at line 1850 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( const Scalar &  alpha) [inherited]

Scale in place: this = alpha*this.

Replace this MultiVector with alpha times this MultiVector. This method will always multiply, even if alpha is zero. That means, for example, that if *this contains NaN entries before calling this method, the NaN entries will remain after this method finishes.

Definition at line 1707 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( Teuchos::ArrayView< const Scalar >  alpha) [inherited]

Scale each column in place: this[j] = alpha[j]*this[j].

Replace each column j of this MultiVector with alpha[j] times the current column j of this MultiVector. This method will always multiply, even if all the entries of alpha are zero. That means, for example, that if *this contains NaN entries before calling this method, the NaN entries will remain after this method finishes.

Definition at line 1737 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A 
) [inherited]

Scale in place: this = alpha * A.

Replace this MultiVector with scaled values of A. This method will always multiply, even if alpha is zero. That means, for example, that if *this contains NaN entries before calling this method, the NaN entries will remain after this method finishes. It is legal for the input A to alias this MultiVector.

Definition at line 1808 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::update ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta 
) [inherited]

Update: this = beta*this + alpha*A.

Update this MultiVector with scaled values of A. If beta is zero, overwrite *this unconditionally, even if it contains NaN entries. It is legal for the input A to alias this MultiVector.

Definition at line 1933 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::update ( const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
const Scalar &  gamma 
) [inherited]

Update: this = gamma*this + alpha*A + beta*B.

Update this MultiVector with scaled values of A and B. If gamma is zero, overwrite *this unconditionally, even if it contains NaN entries. It is legal for the inputs A or B to alias this MultiVector.

Definition at line 1976 of file Tpetra_MultiVector_def.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm1 ( const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &  norms) const [inherited]

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[k] is the one-norm of column k of this MultiVector.

Definition at line 1425 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::norm2 ( const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &  norms) const [inherited]

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 1295 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normInf ( const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &  norms) const [inherited]

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[k] is the infinity-norm of column k of this MultiVector.

Definition at line 1475 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::normWeighted ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  weights,
const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &  norms 
) const [inherited]

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 1351 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::meanValue ( const Teuchos::ArrayView< Scalar > &  means) const [inherited]

Compute mean (average) value 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 1520 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::multiply ( Teuchos::ETransp  transA,
Teuchos::ETransp  transB,
const Scalar &  alpha,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
const Scalar &  beta 
) [inherited]

Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).

If beta is zero, overwrite *this unconditionally, even if it contains NaN entries. This imitates the semantics of analogous BLAS routines like DGEMM.

Definition at line 2724 of file Tpetra_MultiVector_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::elementWiseMultiply ( Scalar  scalarAB,
const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
Scalar  scalarThis 
) [inherited]

Multiply a Vector A elementwise by a MultiVector B.

Compute this = scalarThis * this + scalarAB * B @ A where </tt> denotes element-wise multiplication. In pseudocode, if C denotes *this MultiVector:

 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 2866 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumVectors ( ) const [inline, inherited]

Number of columns in the multivector.

Definition at line 1231 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalLength ( ) const [inherited]

Local number of rows on the calling process.

Definition at line 503 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalLength ( ) const [inherited]

Global number of rows in the multivector.

Definition at line 514 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getStride ( ) const [inherited]

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 525 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isConstantStride ( ) const [inherited]

Whether this multivector has constant stride between columns.

Warning:
This may be different on different processes.

Definition at line 496 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  newMap) [virtual, inherited]

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.

Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3263 of file Tpetra_MultiVector_def.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setCopyOrView ( const Teuchos::DataAccess  copyOrView) [inline, inherited]

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

Warning:
This method is only for expert use. It may change or disappear at any time.

Definition at line 1105 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
Teuchos::DataAccess Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCopyOrView ( ) const [inline, inherited]

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

Warning:
This method is only for expert use. It may change or disappear at any time.

Definition at line 1114 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::assign ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  src) [inherited]

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 3286 of file Tpetra_MultiVector_def.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
template<class T >
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getSubArrayRCP ( ArrayRCP< T >  arr,
size_t  j 
) const [protected, inherited]

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.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::checkSizes ( const SrcDistObject sourceObj) [protected, virtual, inherited]

Whether data redistribution between sourceObj and this object is legal.

This method is called in DistObject::doTransfer() to check whether data redistribution between the two objects is legal.

Implements Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 536 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::constantNumberOfPackets ( ) const [protected, virtual, inherited]

Number of packets to send per LID.

Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 559 of file Tpetra_MultiVector_def.hpp.

virtual void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
size_t  numSameIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  permuteToLIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  permuteFromLIDs 
) [protected, pure virtual, inherited]

Perform copies and permutations that are local to this process.

Parameters:
source[in] On entry, the source object, from which we are distributing. We distribute to the destination object, which is *this object.
numSameIDs[in] The umber of elements that are the same on the source and destination (this) objects. These elements are owned by the same process in both the source and destination objects. No permutation occurs.
numPermuteIDs[in] The number of elements that are locally permuted between the source and destination objects.
permuteToLIDs[in] List of the elements that are permuted. They are listed by their LID in the destination object.
permuteFromLIDs[in] List of the elements that are permuted. They are listed by their LID in the source object.
virtual void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Teuchos::ArrayView< const LocalOrdinal > &  exportLIDs,
Teuchos::Array< Scalar > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets,
Distributor distor 
) [protected, pure virtual, inherited]

Perform any packing or preparation required for communication.

Parameters:
source[in] Source object for the redistribution.
exportLIDs[in] List of the entries (as local IDs in the source object) we will be sending to other images.
exports[out] On exit, the buffer for data to send.
numPacketsPerLID[out] On exit, the implementation of this method must do one of two things: set numPacketsPerLID[i] to contain the number of packets to be exported for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation need not fill numPacketsPerLID.
constantNumPackets[out] On exit, 0 if numPacketsPerLID has variable contents (different size for each LID). If nonzero, then it is expected that the number of packets per LID is constant, and that constantNumPackets is that value.
distor[in] The Distributor object we are using.
virtual void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Teuchos::ArrayView< const LocalOrdinal > &  importLIDs,
const Teuchos::ArrayView< const Scalar > &  imports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t  constantNumPackets,
Distributor distor,
CombineMode  CM 
) [protected, pure virtual, inherited]

Perform any unpacking and combining after communication.

Parameters:
importLIDs[in] List of the entries (as LIDs in the destination object) we received from other images.
imports[in] Buffer containing data we received.
numPacketsPerLID[in] If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i].
constantNumPackets[in] If nonzero, then numPacketsPerLID is constant (same value in all entries) and constantNumPackets is that value. If zero, then numPacketsPerLID[i] is the number of packets imported for importLIDs[i].
distor[in] The Distributor object we are using.
CM[in] The combine mode to use when combining the imported entries with existing entries.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::createViews ( ) const [protected, virtual, inherited]

Hook for creating a const view.

doTransfer() calls this on the source object. By default, it does nothing, but the source object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.

Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3220 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::createViewsNonConst ( KokkosClassic::ReadWriteOption  rwo) [protected, virtual, inherited]

Hook for creating a nonconst view.

doTransfer() calls this on the destination (*this) object. By default, it does nothing, but the destination object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.

Parameters:
rwo[in] Whether to create a write-only or a read-and-write view. For Kokkos Node types where compute buffers live in a separate memory space (e.g., in the device memory of a discrete accelerator like a GPU), a write-only view only requires copying from host memory to the compute buffer, whereas a read-and-write view requires copying both ways (once to read, from the compute buffer to host memory, and once to write, back to the compute buffer).

Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3232 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::releaseViews ( ) const [protected, virtual, inherited]

Hook for releasing views.

doTransfer() calls this on both the source and destination objects, once it no longer needs to access that object's data. By default, this method does nothing. Implementations may use this as a hint to free host memory which is a view of a compute buffer, once the host memory view is no longer needed. Some implementations may prefer to mirror compute buffers in host memory; for these implementations, releaseViews() may do nothing.

Reimplemented from Tpetra::DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3243 of file Tpetra_MultiVector_def.hpp.

void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
) [inherited]

Import data into this object using an Import object ("forward mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Import object if you want to do an Import, else use doExport() with a precomputed Export object.

Parameters:
source[in] The "source" object for redistribution.
importer[in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap().
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
) [inherited]

Import data into this object using an Export object ("reverse mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doImport() that takes a precomputed Import object in that case.

Parameters:
source[in] The "source" object for redistribution.
exporter[in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.)
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
CombineMode  CM 
) [inherited]

Export data into this object using an Export object ("forward mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Export object if you want to do an Export, else use doImport() with a precomputed Import object.

Parameters:
source[in] The "source" object for redistribution.
exporter[in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap().
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
CombineMode  CM 
) [inherited]

Export data into this object using an Import object ("reverse mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doExport() that takes a precomputed Export object in that case.

Parameters:
source[in] The "source" object for redistribution.
importer[in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.)
CM[in] How to combine incoming data with the same global index.
bool Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::isDistributed ( ) const [inherited]

Whether this is a globally distributed object.

For a definition of "globally distributed" (and its opposite, "locally replicated"), see the documentation of Map's isDistributed() method.

virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::getMap ( ) const [inline, virtual, inherited]

The Map describing the parallel distribution of this object.

Note that some Tpetra objects might be distributed using multiple Map objects. For example, CrsMatrix has both a row Map and a column Map. It is up to the subclass to decide which Map to use when invoking the DistObject constructor.

Definition at line 316 of file Tpetra_DistObject_decl.hpp.

void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::print ( std::ostream &  os) const [inherited]

Print this object to the given output stream.

We generally assume that all MPI processes can print to the given stream.

virtual void Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::doTransfer ( const SrcDistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  src,
CombineMode  CM,
size_t  numSameIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  permuteToLIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  permuteFromLIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  remoteLIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  exportLIDs,
Distributor distor,
ReverseOption  revOp 
) [protected, virtual, inherited]

Redistribute data across memory images.

Parameters:
src[in] The source object, to redistribute into the target object, which is *this object.
CM[in] The combine mode that describes how to combine values that map to the same global ID on the same process.
permuteToLIDs[in] See copyAndPermute().
permuteFromLIDs[in] See copyAndPermute().
remoteLIDs[in] List of entries (as local IDs) in the destination object to receive from other processes.
exportLIDs[in] See packAndPrepare().
distor[in/out] The Distributor object that knows how to redistribute data.
revOp[in] Whether to do a forward or reverse mode redistribution.

Friends And Related Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > createMultiVectorFromView ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > &  map,
const Teuchos::ArrayRCP< Scalar > &  view,
const size_t  LDA,
const size_t  numVectors 
) [related]

Nonmember MultiVector constructor with view semantics using user-allocated data.

Warning:
This function is not supported for all Kokkos Node types. Specifically, it is not typically supported for GPU accelerator-based nodes like KokkosClassic::ThrustGPUNode.
Parameters:
map[in] The Map describing the distribution of rows of the multivector.
view[in/out] A pointer to column-major dense matrix data. This will be the multivector's data on the calling process. The multivector will use the pointer directly, without copying.
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.

To Kokkos and Tpetra developers: If you add a new Kokkos Node type that is a host Node type (where memory lives in user space, not in a different space as on a GPU), you will need to add a specialization of Tpetra::details::ViewAccepter for your new Node type.

Definition at line 3843 of file Tpetra_KokkosRefactor_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class DeviceType >
RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > createVectorFromView ( const RCP< const Map< LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< DeviceType > > > &  map,
const ArrayRCP< Scalar > &  view 
) [related]

Non-member function to create a Vector from a specified Map.

Non-member function to create a Vector with view semantics using user-allocated data.

This use case is not supported for all nodes. Specifically, it is not typically supported for accelerator-based nodes like KokkosClassic::ThrustGPUNode.

Definition at line 329 of file Tpetra_KokkosRefactor_Vector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createMultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const size_t  numVectors 
) [related]

Nonmember MultiVector constructor: make a MultiVector from a given Map.

Parameters:
map[in] Map describing the distribution of rows of the resulting MultiVector.
numVectors[in] Number of columns of the resulting MultiVector.

Definition at line 1542 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createMultiVectorFromView ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const Teuchos::ArrayRCP< Scalar > &  view,
const size_t  LDA,
const size_t  numVectors 
) [related]

Nonmember MultiVector constructor with view semantics using user-allocated data.

Warning:
This function is not supported for all Kokkos Node types. Specifically, it is not typically supported for GPU accelerator-based nodes like KokkosClassic::ThrustGPUNode.
Parameters:
map[in] The Map describing the distribution of rows of the multivector.
view[in/out] A pointer to column-major dense matrix data. This will be the multivector's data on the calling process. The multivector will use the pointer directly, without copying.
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.

To Kokkos and Tpetra developers: If you add a new Kokkos Node type that is a host Node type (where memory lives in user space, not in a different space as on a GPU), you will need to add a specialization of Tpetra::details::ViewAccepter for your new Node type.

Definition at line 1589 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createVector ( const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map) [related]

Non-member function to create a Vector from a specified Map.

Definition at line 399 of file Tpetra_Vector_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createVectorFromView ( const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const ArrayRCP< Scalar > &  view 
) [related]

Non-member function to create a Vector with view semantics using user-allocated data.

This use case is not supported for all nodes. Specifically, it is not typically supported for accelerator-based nodes like KokkosClassic::ThrustGPUNode.

Definition at line 411 of file Tpetra_Vector_decl.hpp.

template<class DS , class DL , class DG , class DN , class SS , class SL , class SG , class SN >
void deep_copy ( MultiVector< DS, DL, DG, DN > &  dst,
const MultiVector< SS, SL, SG, SN > &  src 
) [related]

Copy the contents of the MultiVector src into dst.

Precondition:
The two inputs must have the same communicator.
The Map of src must be compatible with the Map of dst.
The two inputs must have the same number of columns.

Copy the contents of the MultiVector src into the MultiVector dst. ("Copy the contents" means the same thing as "deep copy.") The two MultiVectors need not necessarily have the same template parameters, but the assignment of their entries must make sense. Furthermore, their Maps must be compatible, that is, the MultiVectors' local dimensions must be the same on all processes.

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, and throw if they are not compatible.

template<class ST , class LO , class GO , class NT >
MultiVector< ST, LO, GO, NT > createCopy ( const MultiVector< ST, LO, GO, NT > &  src) [related]

Return a deep copy of the MultiVector src.

Regarding Copy or View semantics: The returned MultiVector is always a deep copy of src, and always has view semantics. This is because createCopy returns by value, and therefore it assumes that you want to pass MultiVector objects around by value, not by Teuchos::RCP.

In the Kokkos refactor version of Tpetra, MultiVector always has View semantics. However, the above remarks still apply.


Member Data Documentation

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
KMV Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::lclMV_ [protected, inherited]

The KokkosClassic::MultiVector containing the compute buffer of data.

Definition at line 1140 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
Array<size_t> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::whichVectors_ [protected, inherited]

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 1154 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
bool Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasViewSemantics_ [protected, inherited]

Whether this MultiVector has view semantics.

"View semantics" means that if this MultiVector is on the right side of an operator=, the left side gets a shallow copy, and acquires view semantics. The Kokkos refactor version of MultiVector only ever has view semantics. The "classic" version of MultiVector currently does not have view semantics by default, but this will change.

You can set this for now by calling one of the constructors that accepts a Teuchos::DataAccess enum value.

Definition at line 1167 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
ArrayRCP<Scalar> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ncview_ [mutable, protected, inherited]

Nonconst host view created in createViewsNonConst().

Definition at line 1314 of file Tpetra_MultiVector_decl.hpp.

template<class Scalar = double, class LocalOrdinal = Map<>::local_ordinal_type, class GlobalOrdinal = typename Map<LocalOrdinal>::global_ordinal_type, class Node = typename Map<LocalOrdinal, GlobalOrdinal>::node_type>
ArrayRCP<const Scalar> Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::cview_ [mutable, protected, inherited]

Const host view created in createViews().

Definition at line 1317 of file Tpetra_MultiVector_decl.hpp.

Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::map_ [protected, inherited]

The Map over which this object is distributed.

Definition at line 611 of file Tpetra_DistObject_decl.hpp.

Teuchos::Array<Scalar > Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::imports_ [protected, inherited]

Buffer into which packed data are imported (received from other processes).

Definition at line 615 of file Tpetra_DistObject_decl.hpp.

Teuchos::Array<size_t> Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::numImportPacketsPerLID_ [protected, inherited]

Number of packets to receive for each receive operation.

This array is used in Distributor::doPosts() (and doReversePosts()) when starting the ireceive operation.

This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)

Definition at line 627 of file Tpetra_DistObject_decl.hpp.

Teuchos::Array<Scalar > Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::exports_ [protected, inherited]

Buffer from which packed data are exported (sent to other processes).

Definition at line 630 of file Tpetra_DistObject_decl.hpp.

Teuchos::Array<size_t> Tpetra::DistObject< Scalar , LocalOrdinal, GlobalOrdinal, Node >::numExportPacketsPerLID_ [protected, inherited]

Number of packets to send for each send operation.

This array is used in Distributor::doPosts() (and doReversePosts()) for preparing for the send operation.

This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)

Definition at line 642 of file Tpetra_DistObject_decl.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines