Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Public Member Functions | List of all members
Zoltan2::XpetraMultiVectorAdapter< User > Class Template Reference

An adapter for Xpetra::MultiVector. More...

#include <Zoltan2_XpetraMultiVectorAdapter.hpp>

Inheritance diagram for Zoltan2::XpetraMultiVectorAdapter< User >:
Inheritance graph
[legend]
Collaboration diagram for Zoltan2::XpetraMultiVectorAdapter< User >:
Collaboration graph
[legend]

Public Member Functions

 XpetraMultiVectorAdapter (const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
 Constructor. More...
 
 XpetraMultiVectorAdapter (const RCP< const User > &invector)
 Constructor for case when weights are not being used. More...
 
size_t getLocalNumIDs () const
 Returns the number of objects on this process. More...
 
void getIDsView (const gno_t *&ids) const
 
void getIDsKokkosView (Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
 
int getNumWeightsPerID () const
 Returns the number of weights per object. Number of weights per object should be zero or greater. If zero, then it is assumed that all objects are equally weighted. Default is zero weights per ID. More...
 
void getWeightsView (const scalar_t *&weights, int &stride, int idx) const
 Provide pointer to a weight array with stride. More...
 
void getWeightsKokkos2dView (Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
 
int getNumEntriesPerID () const
 Return the number of vectors. More...
 
void getEntriesView (const scalar_t *&elements, int &stride, int idx=0) const
 Provide a pointer to the elements of the specified vector. More...
 
void getEntriesKokkosView (Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &elements) const
 
template<typename Adapter >
void applyPartitioningSolution (const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
 
template<typename Adapter >
void applyPartitioningSolution (const User &in, RCP< User > &out, const PartitioningSolution< Adapter > &solution) const
 
- Public Member Functions inherited from Zoltan2::VectorAdapter< User >
enum BaseAdapterType adapterType () const override
 Returns the type of adapter. More...
 
virtual void getEntriesKokkosView (Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &elements) const
 Provide a Kokkos view to the elements of the specified vector. More...
 
virtual void getEntriesHostView (typename AdapterWithCoords< User >::CoordsHostView &elements) const
 Provide a Kokkos view (Host side) to the elements of the specified vector. More...
 
virtual void getEntriesDeviceView (typename AdapterWithCoords< User >::CoordsDeviceView &elements) const
 Provide a Kokkos view (Device side) to the elements of the specified vector. More...
 
void generateFiles (const char *fileprefix, const Teuchos::Comm< int > &comm) const
 Write files that can be used as input to Zoltan or Zoltan2 driver Creates chaco-formatted input files for coordinates and weights that can be used as input for Zoltan or Zoltan2 drivers. This routine is SERIAL and can be quite slow. It is meant as a debugging tool only, to allow Zoltan developers to replicate performance that applications are seeing using the applicatios' input. More...
 
int getDimension () const
 
void getCoordinatesView (const scalar_t *&elements, int &stride, int idx=0) const override
 
void getCoordinatesKokkosView (typename AdapterWithCoords< User >::CoordsDeviceView &elements) const override
 
void getCoordinatesHostView (typename AdapterWithCoords< User >::CoordsHostView &elements) const override
 
void getCoordinatesDeviceView (typename AdapterWithCoords< User >::CoordsDeviceView &elements) const override
 
- Public Member Functions inherited from Zoltan2::AdapterWithCoords< User >
virtual void getCoordinatesKokkosView (CoordsDeviceView &elements) const =0
 
virtual void getCoordinatesHostView (CoordsHostView &) const
 
virtual void getCoordinatesDeviceView (CoordsDeviceView &elements) const
 
- Public Member Functions inherited from Zoltan2::BaseAdapter< User >
virtual void getIDsView (const gno_t *&ids) const
 Provide a pointer to this process' identifiers. More...
 
virtual void getIDsKokkosView (ConstIdsDeviceView &ids) const
 Provide a Kokkos view to this process' identifiers. More...
 
virtual void getIDsHostView (ConstIdsHostView &hostIds) const
 Provide a Kokkos view (Host side) to this process' identifiers. More...
 
virtual void getIDsDeviceView (ConstIdsDeviceView &deviceIds) const
 Provide a Kokkos view (Device side) to this process' identifiers. More...
 
virtual void getWeightsKokkosView (Kokkos::View< scalar_t **, device_t > &wgt) const
 Provide kokkos view of weights. More...
 
virtual void getWeightsHostView (WeightsHostView1D &hostWgts, int idx=0) const
 Provide a Kokkos view (Host side) of the weights. More...
 
virtual void getWeightsHostView (WeightsHostView &hostWgts) const
 Provide a Kokkos view (Host side) of the weights. More...
 
virtual void getWeightsDeviceView (WeightsDeviceView1D &deviceWgts, int idx=0) const
 Provide a Kokkos view (Device side) of the weights. More...
 
virtual void getWeightsDeviceView (WeightsDeviceView &deviceWgts) const
 Provide a Kokkos view (Device side) of the weights. More...
 
virtual void getPartsView (const part_t *&inputPart) const
 Provide pointer to an array containing the input part assignment for each ID. The input part information may be used for re-partitioning to reduce data movement, or for mapping parts to processes. Adapters may return NULL for this pointer (the default behavior); if NULL is returned, algorithms will assume the rank. More...
 
virtual void getPartsHostView (Kokkos::View< part_t *, host_t > &inputPart) const
 
virtual void getPartsDeviceView (Kokkos::View< part_t *, device_t > &inputPart) const
 
template<typename Adapter >
void applyPartitioningSolution (const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
 Apply a PartitioningSolution to an input. More...
 
- Public Member Functions inherited from Zoltan2::BaseAdapterRoot
virtual ~BaseAdapterRoot ()=default
 

Additional Inherited Members

- Public Types inherited from Zoltan2::AdapterWithCoords< User >
using scalar_t = typename BaseAdapter< User >::scalar_t
 
using device_t = typename BaseAdapter< User >::node_t::device_type
 
using host_t = typename Kokkos::HostSpace::memory_space
 
using CoordsDeviceView = Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t >
 
using CoordsHostView = typename CoordsDeviceView::HostMirror
 
- Public Types inherited from Zoltan2::BaseAdapter< User >
using lno_t = typename InputTraits< User >::lno_t
 
using gno_t = typename InputTraits< User >::gno_t
 
using scalar_t = typename InputTraits< User >::scalar_t
 
using node_t = typename InputTraits< User >::node_t
 
using part_t = typename InputTraits< User >::part_t
 
using offset_t = typename InputTraits< User >::offset_t
 
using host_t = typename Kokkos::HostSpace::memory_space
 
using device_t = typename node_t::device_type
 
using ConstIdsDeviceView = Kokkos::View< const gno_t *, device_t >
 
using ConstIdsHostView = typename ConstIdsDeviceView::HostMirror
 
using IdsDeviceView = Kokkos::View< gno_t *, device_t >
 
using IdsHostView = typename IdsDeviceView::HostMirror
 
using ConstOffsetsDeviceView = Kokkos::View< const offset_t *, device_t >
 
using ConstOffsetsHostView = typename ConstOffsetsDeviceView::HostMirror
 
using OffsetsDeviceView = Kokkos::View< offset_t *, device_t >
 
using OffsetsHostView = typename OffsetsDeviceView::HostMirror
 
using ConstScalarsDeviceView = Kokkos::View< const scalar_t *, device_t >
 
using ConstScalarsHostView = typename ConstScalarsDeviceView::HostMirror
 
using ScalarsDeviceView = Kokkos::View< scalar_t *, device_t >
 
using ScalarsHostView = typename ScalarsDeviceView::HostMirror
 
using ConstWeightsDeviceView1D = Kokkos::View< const scalar_t *, device_t >
 
using ConstWeightsHostView1D = typename ConstWeightsDeviceView1D::HostMirror
 
using WeightsDeviceView1D = Kokkos::View< scalar_t *, device_t >
 
using WeightsHostView1D = typename WeightsDeviceView1D::HostMirror
 
using ConstWeightsDeviceView = Kokkos::View< const scalar_t **, device_t >
 
using ConstWeightsHostView = typename ConstWeightsDeviceView::HostMirror
 
using WeightsDeviceView = Kokkos::View< scalar_t **, device_t >
 
using WeightsHostView = typename WeightsDeviceView::HostMirror
 
- Protected Member Functions inherited from Zoltan2::BaseAdapter< User >
void generateWeightFileOnly (const char *fileprefix, const Teuchos::Comm< int > &comm) const
 

Detailed Description

template<typename User>
class Zoltan2::XpetraMultiVectorAdapter< User >

An adapter for Xpetra::MultiVector.

The template parameter is the user's input object:

The scalar_t type, representing use data such as vector values, is used by Zoltan2 for weights, coordinates, part sizes and quality metrics. Some User types (like Tpetra::CrsMatrix) have an inherent scalar type, and some (like Tpetra::CrsGraph) do not. For such objects, the scalar type is set by Zoltan2 to float. If you wish to change it to double, set the second template parameter to double.

Definition at line 47 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

Constructor & Destructor Documentation

template<typename User >
Zoltan2::XpetraMultiVectorAdapter< User >::XpetraMultiVectorAdapter ( const RCP< const User > &  invector,
std::vector< const scalar_t * > &  weights,
std::vector< int > &  weightStrides 
)

Constructor.

Parameters
invectorthe user's Xpetra, Tpetra or Epetra MultiVector object
weightsa list of pointers to arrays of weights. The number of weights per multivector element is assumed to be weights.size().
weightStridesa list of strides for the weights. The weight for weight index n for multivector element k should be found at weights[n][weightStrides[n] * k]. If weightStrides.size() is zero, it is assumed all strides are one.

The values pointed to the arguments must remain valid for the lifetime of this Adapter.

Definition at line 201 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User >
Zoltan2::XpetraMultiVectorAdapter< User >::XpetraMultiVectorAdapter ( const RCP< const User > &  invector)

Constructor for case when weights are not being used.

Parameters
invectorthe user's Xpetra, Tpetra or Epetra MultiVector object

Definition at line 234 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

Member Function Documentation

template<typename User>
size_t Zoltan2::XpetraMultiVectorAdapter< User >::getLocalNumIDs ( ) const
inlinevirtual

Returns the number of objects on this process.

Objects may be coordinates, graph vertices, matrix rows, etc. They are the objects to be partitioned, ordered, or colored.

Implements Zoltan2::BaseAdapterRoot.

Definition at line 94 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
void Zoltan2::XpetraMultiVectorAdapter< User >::getIDsView ( const gno_t *&  ids) const
inline

Definition at line 96 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
void Zoltan2::XpetraMultiVectorAdapter< User >::getIDsKokkosView ( Kokkos::View< const gno_t *, typename node_t::device_type > &  ids) const
inline

Definition at line 101 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
int Zoltan2::XpetraMultiVectorAdapter< User >::getNumWeightsPerID ( ) const
inlinevirtual

Returns the number of weights per object. Number of weights per object should be zero or greater. If zero, then it is assumed that all objects are equally weighted. Default is zero weights per ID.

Reimplemented from Zoltan2::BaseAdapterRoot.

Definition at line 131 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
void Zoltan2::XpetraMultiVectorAdapter< User >::getWeightsView ( const scalar_t *&  wgt,
int &  stride,
int  idx 
) const
inlinevirtual

Provide pointer to a weight array with stride.

Parameters
wgton return a pointer to the weights for this idx
strideon return, the value such that the nth weight should be found at wgt[n*stride] .
idxthe weight index, zero or greater This function or getWeightsKokkosView must be implemented in derived adapter if getNumWeightsPerID > 0. This function should not be called if getNumWeightsPerID is zero.

Reimplemented from Zoltan2::BaseAdapter< User >.

Definition at line 133 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
void Zoltan2::XpetraMultiVectorAdapter< User >::getWeightsKokkos2dView ( Kokkos::View< scalar_t **, typename node_t::device_type > &  wgt) const
inline

Definition at line 147 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User>
int Zoltan2::XpetraMultiVectorAdapter< User >::getNumEntriesPerID ( ) const
inlinevirtual

Return the number of vectors.

Implements Zoltan2::VectorAdapter< User >.

Definition at line 169 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User >
void Zoltan2::XpetraMultiVectorAdapter< User >::getEntriesView ( const scalar_t *&  elements,
int &  stride,
int  idx = 0 
) const
virtual

Provide a pointer to the elements of the specified vector.

Parameters
elementswill on return point to the vector values corresponding to the global Ids.
stridethe k'th element is located at elements[stride*k]
idxranges from zero to one less than getNumEntriesPerID(), and represents the vector for which data is being requested.

Reimplemented from Zoltan2::VectorAdapter< User >.

Definition at line 251 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User >
void Zoltan2::XpetraMultiVectorAdapter< User >::getEntriesKokkosView ( Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &  elements) const

Definition at line 293 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User >
template<typename Adapter >
void Zoltan2::XpetraMultiVectorAdapter< User >::applyPartitioningSolution ( const User &  in,
User *&  out,
const PartitioningSolution< Adapter > &  solution 
) const

Definition at line 339 of file Zoltan2_XpetraMultiVectorAdapter.hpp.

template<typename User >
template<typename Adapter >
void Zoltan2::XpetraMultiVectorAdapter< User >::applyPartitioningSolution ( const User &  in,
RCP< User > &  out,
const PartitioningSolution< Adapter > &  solution 
) const

Definition at line 363 of file Zoltan2_XpetraMultiVectorAdapter.hpp.


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