Zoltan2
|
GraphAdapter defines the interface for graph-based user data. More...
#include <Zoltan2_GraphAdapter.hpp>
Public Member Functions | |
enum BaseAdapterType | adapterType () const |
Returns the type of adapter. More... | |
virtual | ~GraphAdapter () |
Destructor. More... | |
GraphAdapter () | |
virtual size_t | getLocalNumVertices () const =0 |
Returns the number of vertices on this process. More... | |
virtual size_t | getLocalNumEdges () const =0 |
Returns the number of edges on this process. More... | |
virtual void | getVertexIDsView (const gno_t *&vertexIds) const =0 |
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 int | getNumWeightsPerVertex () const |
Returns the number (0 or greater) of weights per vertex. More... | |
virtual void | getVertexWeightsView (const scalar_t *&weights, int &stride, int=0) const |
Provide a pointer to the vertex weights, if any. More... | |
virtual bool | useDegreeAsVertexWeight (int) const |
Indicate whether vertex weight with index idx should be the global degree of the vertex. More... | |
virtual int | getNumWeightsPerEdge () const |
Returns the number (0 or greater) of edge weights. 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 |
GraphAdapter defines the interface for graph-based user data.
Adapter objects provide access for Zoltan2 to the user's data. Many built-in adapters are already defined for common data structures, such as Tpetra and Epetra objects and C-language pointers to arrays.
Data types:
scalar_t
vertex and edge weights lno_t
local indices and local counts gno_t
global indices and global counts node_t
is a Kokkos Node typeThe Kokkos node type can be safely ignored.
The template parameter User
is a user-defined data type which, through a traits mechanism, provides the actual data types with which the Zoltan2 library will be compiled. User
may be the actual class or structure used by application to represent a vector, or it may be the helper class BasicUserTypes. See InputTraits for more information.
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 99 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Destructor.
Definition at line 132 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Definition at line 135 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Returns the type of adapter.
Implements Zoltan2::BaseAdapter< User >.
Definition at line 128 of file Zoltan2_GraphAdapter.hpp.
|
pure virtual |
Returns the number of vertices on this process.
Implemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
|
pure virtual |
Returns the number of edges on this process.
Implemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
|
pure virtual |
Sets pointers to this process' graph entries.
vertexIds | will on return a pointer to vertex global Ids |
Implemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
|
pure virtual |
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
offsets | is 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. |
adjIds | on return will point to the array of adjacent vertices for for each vertex. |
Implemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
|
inlinevirtual |
Returns the number (0 or greater) of weights per vertex.
Reimplemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 170 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a pointer to the vertex weights, if any.
weights | is the list of weights of the given index for the vertices returned in getVertexIDsView(). |
stride | The k'th weight is located at weights[stride*k] |
idx | ranges from zero to one less than getNumWeightsPerVertex(). |
Definition at line 178 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Indicate whether vertex weight with index idx should be the global degree of the vertex.
Reimplemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 190 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Returns the number (0 or greater) of edge weights.
Reimplemented in Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >, and Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 197 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a pointer to the edge weights, if any.
weights | is the list of weights of the given index for the edges returned in getEdgeView(). |
stride | The k'th weight is located at weights[stride*k] |
idx | ranges from zero to one less than getNumWeightsPerEdge(). |
Definition at line 205 of file Zoltan2_GraphAdapter.hpp.
|
inline |
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.
coordData | is a pointer to a VectorAdapter with the user's coordinate data. |
Definition at line 222 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Indicate whether coordinate information has been set for this MatrixAdapter.
Definition at line 231 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Obtain the coordinate data registered by the user.
Definition at line 236 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Returns the entity to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_EDGE.
Definition at line 247 of file Zoltan2_GraphAdapter.hpp.
|
inline |
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.
Definition at line 256 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Returns the entity that describes adjacencies between the entities to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_EDGE.
Definition at line 278 of file Zoltan2_GraphAdapter.hpp.
|
inline |
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.
Definition at line 287 of file Zoltan2_GraphAdapter.hpp.
|
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 306 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a pointer to this process' identifiers.
ids | will on return point to the list of the global Ids for this process. |
Reimplemented from Zoltan2::BaseAdapter< User >.
Definition at line 313 of file Zoltan2_GraphAdapter.hpp.
|
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 327 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide pointer to a weight array with stride.
wgt | on return a pointer to the weights for this idx |
stride | on return, the value such that the nth weight should be found at wgt[n*stride] . |
idx | the weight index, zero or greater This function 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 334 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Definition at line 348 of file Zoltan2_GraphAdapter.hpp.