Zoltan2
Loading...
Searching...
No Matches
Zoltan2::TpetraCrsGraphAdapter< User, UserCoord > Class Template Reference

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

#include <Zoltan2_TpetraCrsGraphAdapter.hpp>

Inheritance diagram for Zoltan2::TpetraCrsGraphAdapter< User, UserCoord >:
Collaboration diagram for Zoltan2::TpetraCrsGraphAdapter< User, UserCoord >:

Public Member Functions

 TpetraCrsGraphAdapter (const RCP< const User > &graph, int nVtxWeights=0, int nEdgeWeights=0)
 Constructor for graph with no weights or coordinates.
RCP< const User > getUserGraph () const
 Access to user's graph.
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::TpetraRowGraphAdapter< User, User >
 TpetraRowGraphAdapter (const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
 Constructor for graph with no weights or coordinates.
void setWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to weights for the primary entity type.
void setWeightsDevice (typename Base::ConstWeightsDeviceView1D val, int idx)
 Provide a device view of weights for the primary entity type.
void setWeightsHost (typename Base::ConstWeightsHostView1D val, int idx)
 Provide a host view of weights for the primary entity type.
void setVertexWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to vertex weights.
void setVertexWeightsDevice (typename Base::ConstWeightsDeviceView1D val, int idx)
 Provide a device view to vertex weights.
void setVertexWeightsHost (typename Base::ConstWeightsHostView1D val, int idx)
 Provide a host view to vertex weights.
void setWeightIsDegree (int idx)
 Specify an index for which the weight should be the degree of the entity.
void setVertexWeightIsDegree (int idx)
 Specify an index for which the vertex weight should be the degree of the vertex.
void setEdgeWeights (const scalar_t *val, int stride, int idx)
 Provide a pointer to edge weights.
void setEdgeWeightsDevice (typename Base::ConstWeightsDeviceView1D val, int idx)
 Provide a device view to edge weights.
void setEdgeWeightsHost (typename Base::ConstWeightsHostView1D val, int idx)
 Provide a host view to edge weights.
size_t getLocalNumVertices () const override
 Returns the number of vertices on this process.
void getVertexIDsView (const gno_t *&ids) const override
void getVertexIDsDeviceView (typename Base::ConstIdsDeviceView &ids) const override
 Sets pointers to this process' graph entries.
void getVertexIDsHostView (typename Base::ConstIdsHostView &ids) const override
 Sets pointers to this process' graph entries.
size_t getLocalNumEdges () const override
 Returns the number of edges on this process.
void getEdgesView (const offset_t *&offsets, const gno_t *&adjIds) const override
void getEdgesDeviceView (typename Base::ConstOffsetsDeviceView &offsets, typename Base::ConstIdsDeviceView &adjIds) const override
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
void getEdgesHostView (typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &adjIds) const override
 Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
int getNumWeightsPerVertex () const override
 Returns the number (0 or greater) of weights per vertex.
void getVertexWeightsView (const scalar_t *&weights, int &stride, int idx) const override
void getVertexWeightsDeviceView (typename Base::WeightsDeviceView1D &weights, int idx=0) const override
 Provide a device view of the vertex weights, if any.
void getVertexWeightsHostView (typename Base::WeightsHostView1D &weights, int idx=0) const override
 Provide a host view of the vertex weights, if any.
bool useDegreeAsVertexWeight (int idx) const override
 Indicate whether vertex weight with index idx should be the global degree of the vertex.
int getNumWeightsPerEdge () const override
 Returns the number (0 or greater) of edge weights.
void getEdgeWeightsView (const scalar_t *&weights, int &stride, int idx) const override
void getEdgeWeightsDeviceView (typename Base::WeightsDeviceView1D &weights, int idx=0) const override
 Provide a device view of the edge weights, if any.
void getEdgeWeightsHostView (typename Base::WeightsHostView1D &weights, int idx=0) const override
 Provide a host view of the edge weights, if any.
void applyPartitioningSolution (const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Public Member Functions inherited from Zoltan2::GraphAdapter< User, User >
enum BaseAdapterType adapterType () const override
 Returns the type of adapter.
virtual void getVertexIDsView (const gno_t *&vertexIds) const=0
 Sets pointers to this process' graph entries.
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.
virtual void getVertexWeightsView (const scalar_t *&weights, int &stride, int=0) const
 Provide a pointer to the vertex weights, if any.
virtual void getEdgeWeightsView (const scalar_t *&weights, int &stride, int=0) const
 Provide a pointer to the edge weights, if any.
void setCoordinateInput (VectorAdapter< User > *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_.
bool coordinatesAvailable () const
 Indicate whether coordinate information has been set for this MatrixAdapter.
VectorAdapter< User > * getCoordinateInput () const override
 Obtain the coordinate data registered by the user.
enum GraphEntityType getPrimaryEntityType () const
 Returns the entity to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_EDGE.
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_.
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.
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_.
size_t getLocalNumIDs () const override
 Returns the number of objects on this process.
void getIDsView (const gno_t *&Ids) const override
 Provide a pointer to this process' identifiers.
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.
void getWeightsView (const scalar_t *&wgt, int &stride, int idx=0) const override
 Provide pointer to a weight array with stride.
void getWeightsHostView (typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
void getWeightsDeviceView (typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) 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.
virtual void getIDsHostView (ConstIdsHostView &hostIds) const
 Provide a Kokkos view (Host side) to this process' identifiers.
virtual void getIDsDeviceView (ConstIdsDeviceView &deviceIds) const
 Provide a Kokkos view (Device side) to this process' identifiers.
virtual void getWeightsKokkosView (Kokkos::View< scalar_t **, device_t > &wgt) const
 Provide kokkos view of weights.
virtual void getWeightsHostView (WeightsHostView1D &hostWgts, int idx=0) const
 Provide a Kokkos view (Host side) of the weights.
virtual void getWeightsHostView (WeightsHostView &hostWgts) const
 Provide a Kokkos view (Host side) of the weights.
virtual void getWeightsDeviceView (WeightsDeviceView1D &deviceWgts, int idx=0) const
 Provide a Kokkos view (Device side) of the weights.
virtual void getWeightsDeviceView (WeightsDeviceView &deviceWgts) const
 Provide a Kokkos view (Device side) of the weights.
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.
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.
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::TpetraRowGraphAdapter< User, User >
virtual RCP< User > doMigration (const User &from, size_t numLocalRows, const gno_t *myNewRows) const
Protected Member Functions inherited from Zoltan2::BaseAdapter< User >
void generateWeightFileOnly (const char *fileprefix, const Teuchos::Comm< int > &comm) const
Protected Attributes inherited from Zoltan2::TpetraRowGraphAdapter< User, User >
RCP< const User > graph_
Base::ConstOffsetsHostView offsHost_
Base::ConstIdsHostView adjIdsHost_
Base::ConstIdsDeviceView adjIdsDevice_
Base::ConstOffsetsDeviceView offsDevice_
int nWeightsPerVertex_
ArrayRCP< StridedData< lno_t, scalar_t > > vertexWeights_
Base::WeightsDeviceView vertexWeightsDevice_
Base::VtxDegreeHostView vertexDegreeWeightsHost_
int nWeightsPerEdge_
ArrayRCP< StridedData< lno_t, scalar_t > > edgeWeights_
Base::WeightsDeviceView edgeWeightsDevice_

Detailed Description

template<typename User, typename UserCoord = User>
class Zoltan2::TpetraCrsGraphAdapter< User, UserCoord >

Provides access for Zoltan2 to Tpetra::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.

Definition at line 33 of file Zoltan2_TpetraCrsGraphAdapter.hpp.

Constructor & Destructor Documentation

◆ TpetraCrsGraphAdapter()

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

Constructor for graph with no weights or coordinates.

Parameters
ingraphthe Tpetra::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 82 of file Zoltan2_TpetraCrsGraphAdapter.hpp.

Member Function Documentation

◆ getUserGraph()

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

Access to user's graph.

Definition at line 64 of file Zoltan2_TpetraCrsGraphAdapter.hpp.

◆ applyPartitioningSolution() [1/2]

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

Definition at line 129 of file Zoltan2_TpetraCrsGraphAdapter.hpp.

◆ applyPartitioningSolution() [2/2]

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

Definition at line 139 of file Zoltan2_TpetraCrsGraphAdapter.hpp.


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