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 ()
 Destructor. More...
 
 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 setVertexWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to vertex weights. More...
 
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...
 
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
 Returns the number of vertices on this process. More...
 
void getVertexIDsView (const gno_t *&ids) const
 Sets pointers to this process' graph entries. More...
 
size_t getLocalNumEdges () const
 Returns the number of edges on this process. More...
 
void getEdgesView (const offset_t *&offsets, const gno_t *&adjIds) const
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format. More...
 
int getNumWeightsPerVertex () const
 Returns the number (0 or greater) of weights per vertex. More...
 
void getVertexWeightsView (const scalar_t *&weights, int &stride, int idx) const
 
bool useDegreeAsVertexWeight (int idx) const
 Indicate whether vertex weight with index idx should be the global degree of the vertex. More...
 
int getNumWeightsPerEdge () const
 Returns the number (0 or greater) of edge weights. More...
 
void getEdgeWeightsView (const scalar_t *&weights, int &stride, int idx) 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::GraphAdapter< User, UserCoord >
enum BaseAdapterType adapterType () const
 Returns the type of adapter. More...
 
virtual ~GraphAdapter ()
 Destructor. More...
 
 GraphAdapter ()
 
virtual void getVertexWeightsView (const scalar_t *&weights, int &stride, int=0) const
 Provide a pointer to the vertex weights, if any. More...
 
virtual void getEdgeWeightsView (const scalar_t *&weights, int &stride, int=0) const
 Provide a pointer to the edge weights, if any. More...
 
void setCoordinateInput (VectorAdapter< UserCoord > *coordData)
 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
 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 (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 (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
 Returns the number of objects on this process. More...
 
void getIDsView (const gno_t *&Ids) const
 Provide a pointer to this process' identifiers. More...
 
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 *&wgt, int &stride, int idx=0) const
 Provide pointer to a weight array with stride. More...
 
bool useDegreeAsWeight (int idx) const
 
- Public Member Functions inherited from Zoltan2::BaseAdapter< User >
virtual ~BaseAdapter ()
 Destructor. More...
 
virtual void getIDsKokkosView (Kokkos::View< gno_t * > &) const
 Provide a pointer to this process' identifiers. More...
 
virtual void getWeightsKokkosView (Kokkos::View< scalar_t * > &, int=0) const
 Provide pointer to a weight View. More...
 
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...
 
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 ()
 

Additional Inherited Members

- Public Types inherited from Zoltan2::BaseAdapter< User >
typedef InputTraits< User >::lno_t lno_t
 
typedef InputTraits< User >::gno_t gno_t
 
typedef InputTraits< User >
::scalar_t 
scalar_t
 
typedef InputTraits< User >::part_t part_t
 
typedef InputTraits< User >
::offset_t 
offset_t
 
- 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 84 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

Constructor & Destructor Documentation

template<typename User, typename UserCoord = User>
Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >::~XpetraCrsGraphAdapter ( )
inline

Destructor.

Definition at line 102 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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 292 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 364 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()->getNodeElementList()

Definition at line 375 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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 395 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 411 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()->getNodeElementList()

Then by vertex neighbor:

TheGraph->getLocalRowView(vertexNum, neighborList);

Definition at line 427 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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 187 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 191 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Returns the number of vertices on this process.

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 203 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Sets pointers to this process' graph entries.

Parameters
vertexIdswill on return a pointer to vertex global Ids

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 205 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Returns the number of edges on this process.

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 212 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
inlinevirtual

Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.

Parameters
offsetsis an array of size getLocalNumVertices() + 1. The neighboring vertices for vertexId[i] begin at adjIds[offsets[i]]. The last element of offsets is the size of the adjIds array.
adjIdson return will point to the array of adjacent vertices for for each vertex.

Implements Zoltan2::GraphAdapter< User, UserCoord >.

Definition at line 214 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

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

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

Definition at line 220 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Definition at line 222 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

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

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

Definition at line 238 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

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

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

Definition at line 240 of file Zoltan2_XpetraCrsGraphAdapter.hpp.

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

Definition at line 242 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 448 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 472 of file Zoltan2_XpetraCrsGraphAdapter.hpp.


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