Zoltan2
Loading...
Searching...
No Matches
Zoltan2_GraphAdapter.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef _ZOLTAN2_GRAPHADAPTER_HPP_
15#define _ZOLTAN2_GRAPHADAPTER_HPP_
16
17#include <Zoltan2_Adapter.hpp>
19
20namespace Zoltan2 {
21
25
58template <typename User, typename UserCoord = User>
59class GraphAdapter : public AdapterWithCoordsWrapper<User, UserCoord> {
60private:
64 enum GraphEntityType primaryEntityType_ = GRAPH_VERTEX;
65
69 enum GraphEntityType adjacencyEntityType_ = GRAPH_EDGE;
70
74 VectorAdapter<UserCoord> *coordinateInput_ = nullptr;
75
78 bool haveCoordinateInput_ = false;
79
80public:
81#ifndef DOXYGEN_SHOULD_SKIP_THIS
83 using lno_t = typename InputTraits<User>::lno_t;
84 using gno_t = typename InputTraits<User>::gno_t;
85 using node_t = typename InputTraits<User>::node_t;
87 using user_t = User;
88 using userCoord_t = UserCoord;
89 using base_adapter_t = GraphAdapter<User, UserCoord>;
91 using VtxDegreeHostView = Kokkos::View<bool *, Kokkos::HostSpace>;
92 using device_t = typename node_t::device_type;
93#endif
94
95 enum BaseAdapterType adapterType() const override { return GraphAdapterType; }
96
98 // Methods to be defined in derived classes.
99
102 virtual size_t getLocalNumVertices() const = 0;
103
106 virtual size_t getLocalNumEdges() const = 0;
107
111 virtual void getVertexIDsView(const gno_t *&vertexIds) const = 0;
112
117 virtual void
118 getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &vertexIds) const {
120 }
121
125 virtual void
126 getVertexIDsHostView(typename Base::ConstIdsHostView &vertexIds) const {
128 }
129
139 virtual void getEdgesView(const offset_t *&offsets,
140 const gno_t *&adjIds) const = 0;
141
151 virtual void
152 getEdgesDeviceView(typename Base::ConstOffsetsDeviceView &offsets,
153 typename Base::ConstIdsDeviceView &adjIds) const {
155 }
156
165 virtual void getEdgesHostView(typename Base::ConstOffsetsHostView &offsets,
166 typename Base::ConstIdsHostView &adjIds) const {
168 }
169
172 virtual int getNumWeightsPerVertex() const { return 0; }
173
180 virtual void getVertexWeightsView(const scalar_t *&weights, int &stride,
181 int /* idx */ = 0) const {
182 weights = NULL;
183 stride = 0;
185 }
186
192 virtual void
193 getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
194 int /* idx */ = 0) const {
196 }
197
201 virtual void
202 getVertexWeightsDeviceView(typename Base::WeightsDeviceView &weights) const {
204 }
205
211 virtual void
212 getVertexWeightsHostView(typename Base::WeightsHostView1D &weights,
213 int /* idx */ = 0) const {
215 }
216
220 virtual void
221 getVertexWeightsHostView(typename Base::WeightsHostView &weights) const {
223 }
224
228 virtual bool useDegreeAsVertexWeight(int /* idx */) const { return false; }
229
232 virtual int getNumWeightsPerEdge() const { return 0; }
233
240 virtual void getEdgeWeightsView(const scalar_t *&weights, int &stride,
241 int /* idx */ = 0) const {
242 weights = NULL;
243 stride = 0;
245 }
246
252 virtual void
253 getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
254 int /* idx */ = 0) const {
256 }
257
261 virtual void
262 getEdgeWeightsDeviceView(typename Base::WeightsDeviceView &weights) const {
264 }
265
271 virtual void
272 getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights,
273 int /* idx */ = 0) const {
275 }
276
280 virtual void
281 getEdgeWeightsHostView(typename Base::WeightsHostView &weights) const {
283 }
284
294 coordinateInput_ = coordData;
295 haveCoordinateInput_ = true;
296 }
297
301 bool coordinatesAvailable() const { return haveCoordinateInput_; }
302
307 return coordinateInput_;
308 }
309
311 // Implementations of base-class methods
312
317 return this->primaryEntityType_;
318 }
319
325 void setPrimaryEntityType(const std::string &typestr) {
326 if (typestr == "vertex") {
327 this->primaryEntityType_ = GRAPH_VERTEX;
328 this->adjacencyEntityType_ = GRAPH_EDGE;
329 } else if (typestr == "edge") {
330 this->primaryEntityType_ = GRAPH_EDGE;
331 this->adjacencyEntityType_ = GRAPH_VERTEX;
332 } else {
333 AssertCondition(true, "Invalid GraphEntityType (" + typestr +
334 "). Valid values are 'vertex' and 'edge'");
335 }
336 }
337
343 return this->adjacencyEntityType_;
344 }
345
351 void setAdjacencyEntityType(const std::string &typestr) {
352 if (typestr == "vertex") {
353 this->adjacencyEntityType_ = GRAPH_VERTEX;
354 this->primaryEntityType_ = GRAPH_EDGE;
355 } else if (typestr == "edge") {
356 this->adjacencyEntityType_ = GRAPH_EDGE;
357 this->primaryEntityType_ = GRAPH_VERTEX;
358 } else {
359 AssertCondition(true, "Invalid GraphEntityType (" + typestr +
360 "). Valid values are 'vertex' and 'edge'");
361 }
362 }
363
364 // Functions from the BaseAdapter interface
365 size_t getLocalNumIDs() const override {
367 return getLocalNumVertices();
368 else
369 return getLocalNumEdges();
370 }
371
372 void getIDsView(const gno_t *&Ids) const override {
374 "getIDsView not yet supported for graph edges.");
375
376 getVertexIDsView(Ids);
377 }
378
379 void getIDsDeviceView(typename Base::ConstIdsDeviceView &Ids) const override {
381 "getIDsDeviceView not yet supported for graph edges.");
382
384 }
385
386 void getIDsHostView(typename Base::ConstIdsHostView &Ids) const override {
388 "getIDsHostView not yet supported for graph edges.");
389
391 }
392
393 int getNumWeightsPerID() const override {
395 return getNumWeightsPerVertex();
396 else
397 return getNumWeightsPerEdge();
398 }
399
400 void getWeightsView(const scalar_t *&wgt, int &stride,
401 int idx = 0) const override {
402
404 "getWeightsView not yet supported for graph edges.");
405
406 getVertexWeightsView(wgt, stride, idx);
407 }
408
409 void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts,
410 int idx = 0) const override {
412 "getWeightsHostView not yet supported for graph edges.");
413
414 getVertexWeightsHostView(hostWgts, idx);
415 }
416
417 void getWeightsHostView(typename Base::WeightsHostView &hostWgts) const override {
419 "getWeightsHostView not yet supported for graph edges.");
420
421 getVertexWeightsHostView(hostWgts);
422 }
423
424 void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts,
425 int idx = 0) const override {
427 "getWeightsDeviceView not yet supported for graph edges.");
428
429 getVertexWeightsDeviceView(deviceWgts, idx);
430 }
431
432 void getWeightsDeviceView(typename Base::WeightsDeviceView &deviceWgts) const override {
434 "getWeightsDeviceView not yet supported for graph edges.");
435
436 getVertexWeightsDeviceView(deviceWgts);
437 }
438
439 bool useDegreeAsWeight(int idx) const {
441 "useDegreeAsWeight not yet supported for graph edges.");
442
443 return useDegreeAsVertexWeight(idx);
444 }
445};
446
447} // namespace Zoltan2
448
449#endif
#define Z2_THROW_NOT_IMPLEMENTED
Defines the VectorAdapter interface.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< User >::gno_t gno_t
typename InputTraits< User >::offset_t offset_t
typename node_t::device_type device_t
GraphAdapter defines the interface for graph-based user data.
virtual size_t getLocalNumVertices() const =0
Returns the number of vertices on this process.
virtual void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView &weights) const
Provide a device view of the edge weights, if any.
void setCoordinateInput(VectorAdapter< UserCoord > *coordData) override
Allow user to provide additional data that contains coordinate info associated with the MatrixAdapter...
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const override
virtual void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int=0) const
Provide a device view of the edge 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 getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
virtual void getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int=0) const
Provide a device view of the vertex weights, if any.
virtual bool useDegreeAsVertexWeight(int) const
Indicate whether vertex weight with index idx should be the global degree of the vertex.
void getWeightsDeviceView(typename Base::WeightsDeviceView &deviceWgts) const override
void setPrimaryEntityType(const std::string &typestr)
Sets the primary entity type. Called by algorithm based on parameter value in parameter list from app...
VectorAdapter< UserCoord > * getCoordinateInput() const override
Obtain the coordinate data registered by the user.
void getIDsDeviceView(typename Base::ConstIdsDeviceView &Ids) const override
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.
virtual void getVertexWeightsDeviceView(typename Base::WeightsDeviceView &weights) const
Provide a device view of the vertex weights, if any.
virtual void getVertexWeightsHostView(typename Base::WeightsHostView &weights) const
Provide a host view of the vertex weights, if any.
size_t getLocalNumIDs() const override
Returns the number of objects on this process.
enum GraphEntityType getPrimaryEntityType() const
Returns the entity to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_E...
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
virtual int getNumWeightsPerEdge() const
Returns the number (0 or greater) of edge weights.
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.
void setAdjacencyEntityType(const std::string &typestr)
Sets the adjacency entity type. Called by algorithm based on parameter value in parameter list from a...
void getWeightsHostView(typename Base::WeightsHostView &hostWgts) const override
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 getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &vertexIds) const
Sets pointers to this process' graph entries.
virtual void getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights, int=0) const
Provide a host view of the edge weights, if any.
virtual void getEdgeWeightsHostView(typename Base::WeightsHostView &weights) const
Provide a host view of the edge weights, if any.
virtual void getVertexIDsHostView(typename Base::ConstIdsHostView &vertexIds) const
Sets pointers to this process' graph entries.
virtual void getVertexIDsView(const gno_t *&vertexIds) const =0
Sets pointers to this process' graph entries.
virtual void getVertexWeightsView(const scalar_t *&weights, int &stride, int=0) const
Provide a pointer to the vertex weights, if any.
int getNumWeightsPerID() const override
Returns the number of weights per object. Number of weights per object should be zero or greater....
virtual void getVertexWeightsHostView(typename Base::WeightsHostView1D &weights, int=0) const
Provide a host view of the vertex weights, if any.
bool coordinatesAvailable() const
Indicate whether coordinate information has been set for this MatrixAdapter.
enum GraphEntityType getAdjacencyEntityType() const
Returns the entity that describes adjacencies between the entities to be partitioned,...
void getIDsView(const gno_t *&Ids) const override
Provide a pointer to this process' identifiers.
bool useDegreeAsWeight(int idx) const
virtual int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
enum BaseAdapterType adapterType() const override
Returns the type of adapter.
void getIDsHostView(typename Base::ConstIdsHostView &Ids) const override
virtual size_t getLocalNumEdges() const =0
Returns the number of edges on this process.
VectorAdapter defines the interface for vector input.
Created by mbenlioglu on Aug 31, 2020.
BaseAdapterType
An enum to identify general types of adapters.
@ GraphAdapterType
graph data
GraphEntityType
Enumerated entity type for graphs: Vertices or Edges.
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)
static ArrayRCP< ArrayRCP< zscalar_t > > weights
default_offset_t offset_t
The data type to represent offsets.
default_gno_t gno_t
The ordinal type (e.g., int, long, int64_t) that can represent global counts and identifiers.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices.
default_scalar_t scalar_t
The data type for weights and coordinates.