Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_GraphAdapter.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
50 #ifndef _ZOLTAN2_GRAPHADAPTER_HPP_
51 #define _ZOLTAN2_GRAPHADAPTER_HPP_
52 
53 #include <Zoltan2_Adapter.hpp>
55 
56 namespace Zoltan2 {
57 
61 
94 template <typename User, typename UserCoord = User>
95 class GraphAdapter : public AdapterWithCoordsWrapper<User, UserCoord> {
96 private:
100  enum GraphEntityType primaryEntityType_ = GRAPH_VERTEX;
101 
105  enum GraphEntityType adjacencyEntityType_ = GRAPH_EDGE;
106 
110  VectorAdapter<UserCoord> *coordinateInput_ = nullptr;
111 
114  bool haveCoordinateInput_ = false;
115 
116 public:
117 #ifndef DOXYGEN_SHOULD_SKIP_THIS
118  using scalar_t = typename InputTraits<User>::scalar_t;
119  using lno_t = typename InputTraits<User>::lno_t;
120  using gno_t = typename InputTraits<User>::gno_t;
121  using node_t = typename InputTraits<User>::node_t;
122  using offset_t = typename InputTraits<User>::offset_t;
123  using user_t = User;
124  using userCoord_t = UserCoord;
127  using VtxDegreeHostView = Kokkos::View<bool *, Kokkos::HostSpace>;
128  using device_t = typename node_t::device_type;
129 #endif
130 
131  enum BaseAdapterType adapterType() const override { return GraphAdapterType; }
132 
134  // Methods to be defined in derived classes.
135 
138  virtual size_t getLocalNumVertices() const = 0;
139 
142  virtual size_t getLocalNumEdges() const = 0;
143 
147  virtual void getVertexIDsView(const gno_t *&vertexIds) const = 0;
148 
153  virtual void
154  getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &vertexIds) const {
156  }
157 
161  virtual void
162  getVertexIDsHostView(typename Base::ConstIdsHostView &vertexIds) const {
164  }
165 
175  virtual void getEdgesView(const offset_t *&offsets,
176  const gno_t *&adjIds) const = 0;
177 
187  virtual void
188  getEdgesDeviceView(typename Base::ConstOffsetsDeviceView &offsets,
189  typename Base::ConstIdsDeviceView &adjIds) const {
191  }
192 
201  virtual void getEdgesHostView(typename Base::ConstOffsetsHostView &offsets,
202  typename Base::ConstIdsHostView &adjIds) const {
204  }
205 
208  virtual int getNumWeightsPerVertex() const { return 0; }
209 
216  virtual void getVertexWeightsView(const scalar_t *&weights, int &stride,
217  int /* idx */ = 0) const {
218  weights = NULL;
219  stride = 0;
221  }
222 
228  virtual void
229  getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
230  int /* idx */ = 0) const {
232  }
233 
237  virtual void
238  getVertexWeightsDeviceView(typename Base::WeightsDeviceView &weights) const {
240  }
241 
247  virtual void
248  getVertexWeightsHostView(typename Base::WeightsHostView1D &weights,
249  int /* idx */ = 0) const {
251  }
252 
256  virtual void
257  getVertexWeightsHostView(typename Base::WeightsHostView &weights) const {
259  }
260 
264  virtual bool useDegreeAsVertexWeight(int /* idx */) const { return false; }
265 
268  virtual int getNumWeightsPerEdge() const { return 0; }
269 
276  virtual void getEdgeWeightsView(const scalar_t *&weights, int &stride,
277  int /* idx */ = 0) const {
278  weights = NULL;
279  stride = 0;
281  }
282 
288  virtual void
289  getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
290  int /* idx */ = 0) const {
292  }
293 
297  virtual void
298  getEdgeWeightsDeviceView(typename Base::WeightsDeviceView &weights) const {
300  }
301 
307  virtual void
308  getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights,
309  int /* idx */ = 0) const {
311  }
312 
316  virtual void
317  getEdgeWeightsHostView(typename Base::WeightsHostView &weights) const {
319  }
320 
329  void setCoordinateInput(VectorAdapter<UserCoord> *coordData) override {
330  coordinateInput_ = coordData;
331  haveCoordinateInput_ = true;
332  }
333 
337  bool coordinatesAvailable() const { return haveCoordinateInput_; }
338 
343  return coordinateInput_;
344  }
345 
347  // Implementations of base-class methods
348 
352  inline enum GraphEntityType getPrimaryEntityType() const {
353  return this->primaryEntityType_;
354  }
355 
361  void setPrimaryEntityType(const std::string &typestr) {
362  if (typestr == "vertex") {
363  this->primaryEntityType_ = GRAPH_VERTEX;
364  this->adjacencyEntityType_ = GRAPH_EDGE;
365  } else if (typestr == "edge") {
366  this->primaryEntityType_ = GRAPH_EDGE;
367  this->adjacencyEntityType_ = GRAPH_VERTEX;
368  } else {
369  AssertCondition(true, "Invalid GraphEntityType (" + typestr +
370  "). Valid values are 'vertex' and 'edge'");
371  }
372  }
373 
379  return this->adjacencyEntityType_;
380  }
381 
387  void setAdjacencyEntityType(const std::string &typestr) {
388  if (typestr == "vertex") {
389  this->adjacencyEntityType_ = GRAPH_VERTEX;
390  this->primaryEntityType_ = GRAPH_EDGE;
391  } else if (typestr == "edge") {
392  this->adjacencyEntityType_ = GRAPH_EDGE;
393  this->primaryEntityType_ = GRAPH_VERTEX;
394  } else {
395  AssertCondition(true, "Invalid GraphEntityType (" + typestr +
396  "). Valid values are 'vertex' and 'edge'");
397  }
398  }
399 
400  // Functions from the BaseAdapter interface
401  size_t getLocalNumIDs() const override {
403  return getLocalNumVertices();
404  else
405  return getLocalNumEdges();
406  }
407 
408  void getIDsView(const gno_t *&Ids) const override {
410  "getIDsView not yet supported for graph edges.");
411 
412  getVertexIDsView(Ids);
413  }
414 
415  void getIDsDeviceView(typename Base::ConstIdsDeviceView &Ids) const override {
417  "getIDsDeviceView not yet supported for graph edges.");
418 
420  }
421 
422  void getIDsHostView(typename Base::ConstIdsHostView &Ids) const override {
424  "getIDsHostView not yet supported for graph edges.");
425 
427  }
428 
429  int getNumWeightsPerID() const override {
431  return getNumWeightsPerVertex();
432  else
433  return getNumWeightsPerEdge();
434  }
435 
436  void getWeightsView(const scalar_t *&wgt, int &stride,
437  int idx = 0) const override {
438 
440  "getWeightsView not yet supported for graph edges.");
441 
442  getVertexWeightsView(wgt, stride, idx);
443  }
444 
445  void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts,
446  int idx = 0) const override {
448  "getWeightsHostView not yet supported for graph edges.");
449 
450  getVertexWeightsHostView(hostWgts, idx);
451  }
452 
453  void getWeightsHostView(typename Base::WeightsHostView &hostWgts) const override {
455  "getWeightsHostView not yet supported for graph edges.");
456 
457  getVertexWeightsHostView(hostWgts);
458  }
459 
460  void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts,
461  int idx = 0) const override {
463  "getWeightsDeviceView not yet supported for graph edges.");
464 
465  getVertexWeightsDeviceView(deviceWgts, idx);
466  }
467 
468  void getWeightsDeviceView(typename Base::WeightsDeviceView &deviceWgts) const override {
470  "getWeightsDeviceView not yet supported for graph edges.");
471 
472  getVertexWeightsDeviceView(deviceWgts);
473  }
474 
475  bool useDegreeAsWeight(int idx) const {
477  "useDegreeAsWeight not yet supported for graph edges.");
478 
479  return useDegreeAsVertexWeight(idx);
480  }
481 };
482 
483 } // namespace Zoltan2
484 
485 #endif
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
VectorAdapter< UserCoord > * getCoordinateInput() const override
Obtain the coordinate data registered by the user.
virtual void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int=0) const
Provide a device view of the edge weights, if any.
virtual void getVertexWeightsView(const scalar_t *&weights, int &stride, int=0) const
Provide a pointer to the vertex weights, if any.
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
virtual void getVertexWeightsDeviceView(typename Base::WeightsDeviceView &weights) const
Provide a device view of the vertex weights, if any.
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.
typename InputTraits< User >::scalar_t scalar_t
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)
size_t getLocalNumIDs() const override
Returns the number of objects on this process.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
virtual size_t getLocalNumEdges() const =0
Returns the number of edges on this process.
GraphAdapter defines the interface for graph-based user data.
virtual int getNumWeightsPerEdge() const
Returns the number (0 or greater) of edge weights.
virtual bool useDegreeAsVertexWeight(int) const
Indicate whether vertex weight with index idx should be the global degree of the vertex.
GraphEntityType
Enumerated entity type for graphs: Vertices or Edges.
default_offset_t offset_t
The data type to represent offsets.
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
bool coordinatesAvailable() const
Indicate whether coordinate information has been set for this MatrixAdapter.
virtual void getVertexWeightsHostView(typename Base::WeightsHostView1D &weights, int=0) const
Provide a host view of the vertex weights, if any.
Defines the VectorAdapter interface.
#define Z2_THROW_NOT_IMPLEMENTED
typename node_t::device_type device_t
void getIDsHostView(typename Base::ConstIdsHostView &Ids) const override
virtual void getVertexWeightsHostView(typename Base::WeightsHostView &weights) const
Provide a host view of the vertex weights, if any.
virtual void getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &vertexIds) const
Sets pointers to this process&#39; graph entries.
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const override
Provide pointer to a weight array with stride.
int getNumWeightsPerID() const override
Returns the number of weights per object. Number of weights per object should be zero or greater...
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.
void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
virtual void getVertexIDsHostView(typename Base::ConstIdsHostView &vertexIds) const
Sets pointers to this process&#39; graph entries.
void setCoordinateInput(VectorAdapter< UserCoord > *coordData) override
Allow user to provide additional data that contains coordinate info associated with the MatrixAdapter...
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::gno_t gno_t
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
BaseAdapterType
An enum to identify general types of adapters.
void getWeightsHostView(typename Base::WeightsHostView &hostWgts) const override
virtual void getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int=0) const
Provide a device view of the vertex weights, if any.
void setPrimaryEntityType(const std::string &typestr)
Sets the primary entity type. Called by algorithm based on parameter value in parameter list from app...
void getIDsView(const gno_t *&Ids) const override
Provide a pointer to this process&#39; identifiers.
virtual void getVertexIDsView(const gno_t *&vertexIds) const =0
Sets pointers to this process&#39; graph entries.
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.
virtual int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
typename InputTraits< User >::offset_t offset_t
bool useDegreeAsWeight(int idx) const
virtual void getEdgeWeightsHostView(typename Base::WeightsHostView &weights) const
Provide a host view of the edge weights, if any.
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.
enum GraphEntityType getPrimaryEntityType() const
Returns the entity to be partitioned, ordered, colored, etc. Valid values are GRAPH_VERTEX or GRAPH_E...
virtual size_t getLocalNumVertices() const =0
Returns the number of vertices on this process.
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.
void setAdjacencyEntityType(const std::string &typestr)
Sets the adjacency entity type. Called by algorithm based on parameter value in parameter list from a...
void getIDsDeviceView(typename Base::ConstIdsDeviceView &Ids) const override
virtual void getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights, int=0) const
Provide a host view of the edge weights, if any.
typename InputTraits< User >::lno_t lno_t
default_scalar_t scalar_t
The data type for weights and coordinates.
virtual void getEdgeWeightsView(const scalar_t *&weights, int &stride, int=0) const
Provide a pointer to the edge weights, if any.
virtual void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView &weights) const
Provide a device view of the edge weights, if any.
enum BaseAdapterType adapterType() const override
Returns the type of adapter.
void getWeightsDeviceView(typename Base::WeightsDeviceView &deviceWgts) const override