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

Provides access for Zoltan2 to Xpetra::CrsGraph data. More...

#include <Zoltan2_XpetraCrsGraphAdapter.hpp>

Inheritance diagram for Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >:
Inheritance graph
[legend]
Collaboration diagram for Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >:
Collaboration graph
[legend]

Public Member Functions

 XpetraCrsGraphAdapter (const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
 Constructor for graph with no weights or coordinates. More...
 
void setWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to weights for the primary entity type. More...
 
void setWeightsDevice (typename Base::ConstWeightsDeviceView &val, int idx)
 
void setWeightsHost (typename Base::ConstWeightsHostView &val, int idx)
 
void setVertexWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to vertex weights. More...
 
void setVertexWeightsDevice (typename Base::ConstWeightsDeviceView &val, int idx)
 
void setVertexWeightsHost (typename Base::ConstWeightsHostView &val, int idx)
 
void setWeightIsDegree (int idx)
 Specify an index for which the weight should be the degree of the entity. More...
 
void setVertexWeightIsDegree (int idx)
 Specify an index for which the vertex weight should be the degree of the vertex. More...
 
void setEdgeWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to edge weights. More...
 
void setEdgeWeightsDevice (typename Base::ConstWeightsDeviceView &val, int idx)
 
void setEdgeWeightsHost (typename Base::ConstWeightsHostView &val, int idx)
 
RCP< const xgraph_tgetXpetraGraph () const
 Access to Xpetra-wrapped user's graph. More...
 
RCP< const User > getUserGraph () const
 Access to user's graph. More...
 
size_t getLocalNumVertices () const override
 Returns the number of vertices on this process. More...
 
void getVertexIDsView (const gno_t *&ids) const override
 
size_t getLocalNumEdges () const override
 Returns the number of edges on this process. More...
 
void getEdgesView (const offset_t *&offsets, const gno_t *&adjIds) const override
 
int getNumWeightsPerVertex () const override
 Returns the number (0 or greater) of weights per vertex. More...
 
void getVertexWeightsView (const scalar_t *&weights, int &stride, int idx) const override
 Provide a pointer to the vertex weights, if any. More...
 
bool useDegreeAsVertexWeight (int idx) const override
 Indicate whether vertex weight with index idx should be the global degree of the vertex. More...
 
int getNumWeightsPerEdge () const override
 Returns the number (0 or greater) of edge weights. More...
 
void getEdgeWeightsView (const scalar_t *&weights, int &stride, int idx) const override
 Provide a pointer to the edge weights, if any. More...
 
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::GraphAdapter< User, UserCoord >
enum BaseAdapterType adapterType () const override
 Returns the type of adapter. More...
 
virtual void getVertexIDsView (const gno_t *&vertexIds) const =0
 Sets pointers to this process' graph entries. More...
 
virtual void getVertexIDsDeviceView (typename Base::ConstIdsDeviceView &vertexIds) const
 Sets pointers to this process' graph entries. More...
 
virtual void getVertexIDsHostView (typename Base::ConstIdsHostView &vertexIds) const
 Sets pointers to this process' graph entries. More...
 
virtual void getEdgesView (const offset_t *&offsets, const gno_t *&adjIds) const =0
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format. More...
 
virtual void getEdgesDeviceView (typename Base::ConstOffsetsDeviceView &offsets, typename Base::ConstIdsDeviceView &adjIds) const
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format. More...
 
virtual void getEdgesHostView (typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &adjIds) const
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format. More...
 
virtual void getVertexWeightsDeviceView (typename Base::WeightsDeviceView1D &weights, int=0) const
 Provide a device view of the vertex weights, if any. More...
 
virtual void getVertexWeightsDeviceView (typename Base::WeightsDeviceView &weights) const
 Provide a device view of the vertex weights, if any. More...
 
virtual void getVertexWeightsHostView (typename Base::WeightsHostView1D &weights, int=0) const
 Provide a host view of the vertex weights, if any. More...
 
virtual void getVertexWeightsHostView (typename Base::WeightsHostView &weights) const
 Provide a host view of the vertex weights, if any. More...
 
virtual void getEdgeWeightsDeviceView (typename Base::WeightsDeviceView1D &weights, int=0) const
 Provide a device view of the edge weights, if any. More...
 
virtual void getEdgeWeightsDeviceView (typename Base::WeightsDeviceView &weights) const
 Provide a device view of the edge weights, if any. More...
 
virtual void getEdgeWeightsHostView (typename Base::WeightsHostView1D &weights, int=0) const
 Provide a host view of the edge weights, if any. More...
 
virtual void getEdgeWeightsHostView (typename Base::WeightsHostView &weights) const
 Provide a host view of the edge weights, if any. More...
 
void setCoordinateInput (VectorAdapter< UserCoord > *coordData) override
 Allow user to provide additional data that contains coordinate info associated with the MatrixAdapter's primaryEntityType_. Associated data must have the same parallel distribution and ordering of entries as the primaryEntityType_. More...
 
bool coordinatesAvailable () const
 Indicate whether coordinate information has been set for this MatrixAdapter. More...
 
VectorAdapter< UserCoord > * getCoordinateInput () const override
 Obtain the coordinate data registered by the user. More...
 
enum GraphEntityType getPrimaryEntityType () const
 Returns the entity to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_EDGE. More...
 
void setPrimaryEntityType (const std::string &typestr)
 Sets the primary entity type. Called by algorithm based on parameter value in parameter list from application. Also sets to adjacencyEntityType_ to something reasonable: opposite of primaryEntityType_. More...
 
enum GraphEntityType getAdjacencyEntityType () const
 Returns the entity that describes adjacencies between the entities to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_EDGE. More...
 
void setAdjacencyEntityType (const std::string &typestr)
 Sets the adjacency entity type. Called by algorithm based on parameter value in parameter list from application. Also sets to primaryEntityType_ to something reasonable: opposite of adjacencyEntityType_. More...
 
size_t getLocalNumIDs () const override
 Returns the number of objects on this process. More...
 
void getIDsView (const gno_t *&Ids) const override
 Provide a pointer to this process' identifiers. More...
 
void getIDsDeviceView (typename Base::ConstIdsDeviceView &Ids) const override
 
void getIDsHostView (typename Base::ConstIdsHostView &Ids) const override
 
int getNumWeightsPerID () const override
 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 *&wgt, int &stride, int idx=0) const override
 Provide pointer to a weight array with stride. More...
 
void getWeightsHostView (typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
 
void getWeightsHostView (typename Base::WeightsHostView &hostWgts) const override
 
void getWeightsDeviceView (typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
 
void getWeightsDeviceView (typename Base::WeightsDeviceView &deviceWgts) const override
 
bool useDegreeAsWeight (int idx) const
 
- Public Member Functions inherited from Zoltan2::BaseAdapter< User >
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::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, typename UserCoord = User>
class Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >

Provides access for Zoltan2 to Xpetra::CrsGraph data.

Todo:

test for memory alloc failure when we resize a vector

we assume FillComplete has been called. Should we support objects that are not FillCompleted.

The template parameter is the user's input object:

The scalar_t type, representing use data such as matrix 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 48 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

Constructor & Destructor Documentation

template<typename User , typename UserCoord >
Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::XpetraCrsGraphAdapter ( const RCP< const User > &  ingraph,
int  nVtxWeights = 0,
int  nEdgeWeights = 0 
)

Constructor for graph with no weights or coordinates.

Parameters
ingraphthe Epetra_CrsGraph, Tpetra::CrsGraph or Xpetra::CrsGraph
numVtxWeightsthe number of weights per vertex (default = 0)
numEdgeWeightsthe number of weights per edge (default = 0)

Most adapters do not have RCPs in their interface. This one does because the user is obviously a Trilinos user.

Definition at line 255 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

Member Function Documentation

template<typename User , typename UserCoord >
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setWeights ( const scalar_t val,
int  stride,
int  idx 
)

Provide a pointer to weights for the primary entity type.

Parameters
valA pointer to the weights for index idx.
strideA stride for the val array. If is k, then val[n * k] is the weight for the n th entity for index .
idxA number from 0 to one less than weight idx specified in the constructor.

The order of the weights should match the order that entities appear in the input data structure.

Definition at line 327 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setWeightsDevice ( typename Base::ConstWeightsDeviceView &  val,
int  idx 
)
inline

Definition at line 91 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setWeightsHost ( typename Base::ConstWeightsHostView &  val,
int  idx 
)
inline

Definition at line 92 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User , typename UserCoord >
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setVertexWeights ( const scalar_t val,
int  stride,
int  idx 
)

Provide a pointer to vertex weights.

Parameters
valA pointer to the weights for index idx.
strideA stride for the val array. If is k, then val[n * k] is the weight for the n th vertex for index .
idxA number from 0 to one less than number of vertex weights specified in the constructor.

The order of the vertex weights should match the order that vertices appear in the input data structure.

TheGraph->getRowMap()->getLocalElementList()

Definition at line 338 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setVertexWeightsDevice ( typename Base::ConstWeightsDeviceView &  val,
int  idx 
)
template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setVertexWeightsHost ( typename Base::ConstWeightsHostView &  val,
int  idx 
)
template<typename User , typename UserCoord >
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setWeightIsDegree ( int  idx)

Specify an index for which the weight should be the degree of the entity.

Parameters
idxZoltan2 will use the entity's degree as the entity weight for index idx.

Definition at line 358 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User , typename UserCoord >
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setVertexWeightIsDegree ( int  idx)

Specify an index for which the vertex weight should be the degree of the vertex.

Parameters
idxZoltan2 will use the vertex's degree as the vertex weight for index idx.

Definition at line 374 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User , typename UserCoord >
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setEdgeWeights ( const scalar_t val,
int  stride,
int  idx 
)

Provide a pointer to edge weights.

Parameters
valA pointer to the weights for index idx.
strideA stride for the val array. If is k, then val[n * k] is the weight for the n th edge for index .
dimA number from 0 to one less than the number of edge weights specified in the constructor.

The order of the edge weights should follow the order that the the vertices and edges appear in the input data structure.

By vertex:

TheGraph->getRowMap()->getLocalElementList()

Then by vertex neighbor:

TheGraph->getLocalRowView(vertexNum, neighborList);

Definition at line 390 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setEdgeWeightsDevice ( typename Base::ConstWeightsDeviceView &  val,
int  idx 
)
template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::setEdgeWeightsHost ( typename Base::ConstWeightsHostView &  val,
int  idx 
)
template<typename User, typename UserCoord = User>
RCP<const xgraph_t> Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getXpetraGraph ( ) const
inline

Access to Xpetra-wrapped user's graph.

Definition at line 155 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
RCP<const User> Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getUserGraph ( ) const
inline

Access to user's graph.

Definition at line 159 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
size_t Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getLocalNumVertices ( ) const
inlineoverridevirtual

Returns the number of vertices on this process.

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 167 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getVertexIDsView ( const gno_t *&  ids) const
inlineoverride

Definition at line 169 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
size_t Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getLocalNumEdges ( ) const
inlineoverridevirtual

Returns the number of edges on this process.

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 176 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getEdgesView ( const offset_t *&  offsets,
const gno_t *&  adjIds 
) const
inlineoverride

Definition at line 178 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
int Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getNumWeightsPerVertex ( ) const
inlineoverridevirtual

Returns the number (0 or greater) of weights per vertex.

Reimplemented from Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 184 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getVertexWeightsView ( const scalar_t *&  weights,
int &  stride,
int   
) const
inlineoverridevirtual

Provide a pointer to the vertex weights, if any.

Parameters
weightsis the list of weights of the given index for the vertices returned in getVertexIDsView().
strideThe k'th weight is located at weights[stride*k]
idxranges from zero to one less than getNumWeightsPerVertex().

Reimplemented from Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 186 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
bool Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::useDegreeAsVertexWeight ( int  ) const
inlineoverridevirtual

Indicate whether vertex weight with index idx should be the global degree of the vertex.

Reimplemented from Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 202 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
int Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getNumWeightsPerEdge ( ) const
inlineoverridevirtual

Returns the number (0 or greater) of edge weights.

Reimplemented from Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 204 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

template<typename User, typename UserCoord = User>
void Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::getEdgeWeightsView ( const scalar_t *&  weights,
int &  stride,
int   
) const
inlineoverridevirtual

Provide a pointer to the edge weights, if any.

Parameters
weightsis the list of weights of the given index for the edges returned in getEdgeView().
strideThe k'th weight is located at weights[stride*k]
idxranges from zero to one less than getNumWeightsPerEdge().

Reimplemented from Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 206 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Definition at line 411 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Definition at line 435 of file Zoltan2_XpetraCrsGraphAdapter.hpp.


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