Zoltan2
|
GraphAdapter defines the interface for graph-based user data. More...
#include <Zoltan2_GraphAdapter.hpp>
Public Member Functions | |
enum BaseAdapterType | adapterType () const override |
Returns the type of adapter. More... | |
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 | 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 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 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 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... | |
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 |
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 59 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
Returns the type of adapter.
Implements Zoltan2::BaseAdapter< User >.
Definition at line 95 of file Zoltan2_GraphAdapter.hpp.
|
pure virtual |
Returns the number of vertices on this process.
Implemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
|
pure virtual |
Returns the number of edges on this process.
Implemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
|
pure virtual |
Sets pointers to this process' graph entries.
vertexIds | will on return a pointer to vertex global Ids |
Implemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
|
inlinevirtual |
Sets pointers to this process' graph entries.
vertexIds | will on return a device Kokkos::View with vertex global Ids |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 118 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Sets pointers to this process' graph entries.
vertexIds | will on return a host Kokkos::View with vertex global Ids |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 126 of file Zoltan2_GraphAdapter.hpp.
|
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::TpetraRowGraphAdapter< User, UserCoord >.
|
inlinevirtual |
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
offsets | is device Kokkos::View 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 | Device Kokkos::View of adjacent vertices for for each vertex. |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 152 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
offsets | is host Kokkos::View 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 | Host Kokkos::View of adjacent vertices for for each vertex. |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 165 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Returns the number (0 or greater) of weights per vertex.
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
Definition at line 172 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(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
Definition at line 180 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a device view of the vertex weights, if any.
weights | is the list of weights of the given index for the vertices returned in getVertexIDsView(). |
idx | ranges from zero to one less than getNumWeightsPerVertex(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 193 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a device view of the vertex weights, if any.
weights | is the view of all the weights for the vertices returned in getVertexIDsView(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 202 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a host view of the vertex weights, if any.
weights | is the list of weights of the given index for the vertices returned in getVertexIDsView(). |
idx | ranges from zero to one less than getNumWeightsPerVertex(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 212 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a host view of the vertex weights, if any.
weights | is the list of all the weights for the vertices returned in getVertexIDsView() |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 221 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Indicate whether vertex weight with index idx should be the global degree of the vertex.
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
Definition at line 228 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Returns the number (0 or greater) of edge weights.
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
Definition at line 232 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(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >, and Zoltan2::XpetraCrsGraphAdapter< User, UserCoord >.
Definition at line 240 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a device view of the edge weights, if any.
weights | is the list of weights of the given index for the edges returned in getEdgeView(). |
idx | ranges from zero to one less than getNumWeightsPerEdge(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 253 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a device view of the edge weights, if any.
weights | is the list of weights for the edges returned in getEdgeView(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 262 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a host view of the edge weights, if any.
weights | is the list of weights of the given index for the edges returned in getEdgeView(). |
idx | ranges from zero to one less than getNumWeightsPerEdge(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 272 of file Zoltan2_GraphAdapter.hpp.
|
inlinevirtual |
Provide a host view of the edge weights, if any.
weights | is the list of weights for the edges returned in getEdgeView(). |
Reimplemented in Zoltan2::TpetraRowGraphAdapter< User, UserCoord >.
Definition at line 281 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
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. |
Implements Zoltan2::AdapterWithCoordsWrapper< User, UserCoord >.
Definition at line 293 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Indicate whether coordinate information has been set for this MatrixAdapter.
Definition at line 301 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
Obtain the coordinate data registered by the user.
Implements Zoltan2::AdapterWithCoordsWrapper< User, UserCoord >.
Definition at line 306 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 316 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 325 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 342 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 351 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
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 365 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
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 372 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 379 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 386 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
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 393 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverridevirtual |
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 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 400 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 409 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 417 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 424 of file Zoltan2_GraphAdapter.hpp.
|
inlineoverride |
Definition at line 432 of file Zoltan2_GraphAdapter.hpp.
|
inline |
Definition at line 439 of file Zoltan2_GraphAdapter.hpp.