Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Amesos2::MultiVecAdapter< Epetra_MultiVector > Class Template Reference

Amesos2 adapter for the Epetra_MultiVector class. More...

#include <Amesos2_EpetraMultiVecAdapter_decl.hpp>

Public Types

typedef double scalar_t
 
typedef int local_ordinal_t
 
typedef
Tpetra::Map::global_ordinal_type 
global_ordinal_t
 
typedef size_t global_size_t
 
typedef Tpetra::Map::node_type node_t
 
typedef Epetra_MultiVector multivec_t
 

Public Member Functions

bool isLocallyIndexed () const
 Checks whether this multi-vector is local to the calling node.
 
bool isGloballyIndexed () const
 
Teuchos::RCP< const
Tpetra::Map< local_ordinal_t,
global_ordinal_t, node_t > > 
getMap () const
 Get a Tpetra::Map that describes this MultiVector. More...
 
const Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 Returns the Teuchos::Comm object associated with this multi-vector.
 
size_t getLocalLength () const
 Get the length of vectors local to the calling node.
 
size_t getLocalNumVectors () const
 Get the number of vectors on this node.
 
global_size_t getGlobalLength () const
 Get the length of vectors in the global space.
 
size_t getGlobalNumVectors () const
 Get the number of global vectors.
 
size_t getStride () const
 Return the stride between vectors on this node.
 
bool isConstantStride () const
 Return true if this MV has constant stride between vectors on this node.
 
Teuchos::RCP< const
Tpetra::Vector< scalar_t,
local_ordinal_t,
global_ordinal_t, node_t > > 
getVector (size_t j) const
 Const vector access.
 
Teuchos::RCP< Tpetra::Vector
< scalar_t, local_ordinal_t,
global_ordinal_t, node_t > > 
getVectorNonConst (size_t j)
 Nonconst vector access. More...
 
double * getMVPointer_impl () const
 Return pointer to vector when number of vectors == 1 and num procs == 1.
 
void get1dCopy (const Teuchos::ArrayView< scalar_t > &A, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > distribution_map, EDistribution distribution) const
 Copies the multi-vector's data into the user-provided vector. More...
 
template<typename KV >
void get1dCopy_kokkos_view (KV &A, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > distribution_map, EDistribution distribution) const
 
void get1dCopy_kokkos_view_host (Kokkos::View< scalar_t **, Kokkos::LayoutLeft, Kokkos::HostSpace > &new_data, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > distribution_map, EDistribution) const
 
Teuchos::ArrayRCP< scalar_t > get1dViewNonConst (bool local=false)
 Extracts a 1 dimensional view of this multi-vector's data. More...
 
void put1dData (const Teuchos::ArrayView< const scalar_t > &new_data, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > source_map, EDistribution distribution)
 Export newVals into the global MultiVector space. More...
 
template<typename KV >
void put1dData_kokkos_view (KV &new_data, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > source_map, EDistribution distribution)
 
void put1dData_kokkos_view_host (Kokkos::View< scalar_t **, Kokkos::LayoutLeft, Kokkos::HostSpace > &new_data, size_t lda, Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > > source_map, EDistribution distribution)
 
std::string description () const
 Get a short description of this adapter class.
 
void describe (Teuchos::FancyOStream &os, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print a description of this adapter to the Fancy Output Stream.
 

Static Public Attributes

static const char * name = "Amesos2 adapter for Epetra_MultiVector"
 

Protected Member Functions

 MultiVecAdapter (const MultiVecAdapter< multivec_t > &adapter)
 Copy constructor.
 
 MultiVecAdapter (const Teuchos::RCP< multivec_t > &m)
 Initialize an adapter from a multi-vector RCP. More...
 

Private Attributes

Teuchos::RCP< multivec_t > mv_
 The multi-vector this adapter wraps.
 
Teuchos::RCP< Epetra_Import > importer_
 
Teuchos::RCP< Epetra_Export > exporter_
 
Teuchos::RCP< const
Epetra_BlockMap > 
mv_map_
 

Friends

Teuchos::RCP< MultiVecAdapter
< multivec_t > > 
createMultiVecAdapter (Teuchos::RCP< multivec_t >)
 
Teuchos::RCP< const
MultiVecAdapter< multivec_t > > 
createConstMultiVecAdapter (Teuchos::RCP< const multivec_t >)
 

Detailed Description

template<>
class Amesos2::MultiVecAdapter< Epetra_MultiVector >

Amesos2 adapter for the Epetra_MultiVector class.

Constructor & Destructor Documentation

Amesos2::MultiVecAdapter< Epetra_MultiVector >::MultiVecAdapter ( const Teuchos::RCP< multivec_t > &  m)
protected

Initialize an adapter from a multi-vector RCP.

Parameters
mAn RCP pointing to the multi-vector which is to be wrapped.

Member Function Documentation

Teuchos::RCP< const Tpetra::Map< MultiVecAdapter< Epetra_MultiVector >::local_ordinal_t, MultiVecAdapter< Epetra_MultiVector >::global_ordinal_t, MultiVecAdapter< Epetra_MultiVector >::node_t > > Amesos2::MultiVecAdapter< Epetra_MultiVector >::getMap ( ) const

Get a Tpetra::Map that describes this MultiVector.

Not part of the MultiVecAdapter interface, but useful for other adaptations.

Teuchos::RCP< Tpetra::Vector< MultiVecAdapter< Epetra_MultiVector >::scalar_t, MultiVecAdapter< Epetra_MultiVector >::local_ordinal_t, MultiVecAdapter< Epetra_MultiVector >::global_ordinal_t, MultiVecAdapter< Epetra_MultiVector >::node_t > > Amesos2::MultiVecAdapter< Epetra_MultiVector >::getVectorNonConst ( size_t  j)

Nonconst vector access.

Note
Vectors returned hold a copy of the data in the multi-vector. So any changes to the returned vector will not be represented in the underlying multi-vector.
void Amesos2::MultiVecAdapter< Epetra_MultiVector >::get1dCopy ( const Teuchos::ArrayView< scalar_t > &  A,
size_t  lda,
Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > >  distribution_map,
EDistribution  distribution 
) const

Copies the multi-vector's data into the user-provided vector.

Each multi-vector is lda apart in memory.

References Amesos2::Util::tpetra_map_to_epetra_map().

Teuchos::ArrayRCP< MultiVecAdapter< Epetra_MultiVector >::scalar_t > Amesos2::MultiVecAdapter< Epetra_MultiVector >::get1dViewNonConst ( bool  local = false)

Extracts a 1 dimensional view of this multi-vector's data.

Guarantees that the view returned will reside in contiguous storage.

Warning
It is recommended to use the get1dCopy function, from a data-hiding perspective. Use if you know what you are doing.
Parameters
localif true , each node will get a view of the vectors it is in possession of. The default, false , will give each calling node a view of the global multi-vector.
Note
This function is not declared const as it normally would be, since it must modify local copies of the vector data before returning the result.
void Amesos2::MultiVecAdapter< Epetra_MultiVector >::put1dData ( const Teuchos::ArrayView< const scalar_t > &  new_data,
size_t  lda,
Teuchos::Ptr< const Tpetra::Map< local_ordinal_t, global_ordinal_t, node_t > >  source_map,
EDistribution  distribution 
)

Export newVals into the global MultiVector space.

Note
we assume the vectors in newVals have the same leading dimension as those in this
Template Parameters
Value_tThe type of the data values that are being put into mv_
Parameters
newValsThe values to be exported into the global space.

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