Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_XpetraCrsGraphAdapter.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_XPETRACRSGRAPHADAPTER_HPP_
51 #define _ZOLTAN2_XPETRACRSGRAPHADAPTER_HPP_
52 
53 #include <Zoltan2_GraphAdapter.hpp>
54 #include <Zoltan2_StridedData.hpp>
55 #include <Zoltan2_XpetraTraits.hpp>
57 
58 #include <Xpetra_CrsGraph.hpp>
59 
60 namespace Zoltan2 {
61 
83 template <typename User, typename UserCoord=User>
84  class XpetraCrsGraphAdapter : public GraphAdapter<User,UserCoord> {
85 
86 public:
87 
88 #ifndef DOXYGEN_SHOULD_SKIP_THIS
89  typedef typename InputTraits<User>::scalar_t scalar_t;
90  typedef typename InputTraits<User>::offset_t offset_t;
91  typedef typename InputTraits<User>::lno_t lno_t;
92  typedef typename InputTraits<User>::gno_t gno_t;
93  typedef typename InputTraits<User>::part_t part_t;
94  typedef typename InputTraits<User>::node_t node_t;
95  typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
96  typedef User user_t;
97  typedef UserCoord userCoord_t;
98 #endif
99 
103 
113  XpetraCrsGraphAdapter(const RCP<const User> &ingraph,
114  int nVtxWeights=0, int nEdgeWeights=0);
115 
128  void setWeights(const scalar_t *val, int stride, int idx);
129 
145  void setVertexWeights(const scalar_t *val, int stride, int idx);
146 
152  void setWeightIsDegree(int idx);
153 
159  void setVertexWeightIsDegree(int idx);
160 
183  void setEdgeWeights(const scalar_t *val, int stride, int idx);
184 
187  RCP<const xgraph_t> getXpetraGraph() const { return graph_; }
188 
191  RCP<const User> getUserGraph() const { return ingraph_; }
192 
194  // The Adapter interface.
196 
198  // The GraphAdapter interface.
200 
201  // TODO: Assuming rows == objects;
202  // TODO: Need to add option for columns or nonzeros?
203  size_t getLocalNumVertices() const { return graph_->getNodeNumRows(); }
204 
205  void getVertexIDsView(const gno_t *&ids) const
206  {
207  ids = NULL;
208  if (getLocalNumVertices())
209  ids = graph_->getRowMap()->getNodeElementList().getRawPtr();
210  }
211 
212  size_t getLocalNumEdges() const { return graph_->getNodeNumEntries(); }
213 
214  void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
215  {
216  offsets = offs_.getRawPtr();
217  adjIds = (getLocalNumEdges() ? adjids_.getRawPtr() : NULL);
218  }
219 
220  int getNumWeightsPerVertex() const { return nWeightsPerVertex_;}
221 
222  void getVertexWeightsView(const scalar_t *&weights, int &stride,
223  int idx) const
224  {
225  if(idx<0 || idx >= nWeightsPerVertex_)
226  {
227  std::ostringstream emsg;
228  emsg << __FILE__ << ":" << __LINE__
229  << " Invalid vertex weight index " << idx << std::endl;
230  throw std::runtime_error(emsg.str());
231  }
232 
233 
234  size_t length;
235  vertexWeights_[idx].getStridedList(length, weights, stride);
236  }
237 
238  bool useDegreeAsVertexWeight(int idx) const {return vertexDegreeWeight_[idx];}
239 
240  int getNumWeightsPerEdge() const { return nWeightsPerEdge_;}
241 
242  void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
243  {
244  if(idx<0 || idx >= nWeightsPerEdge_)
245  {
246  std::ostringstream emsg;
247  emsg << __FILE__ << ":" << __LINE__
248  << " Invalid edge weight index " << idx << std::endl;
249  throw std::runtime_error(emsg.str());
250  }
251 
252 
253  size_t length;
254  edgeWeights_[idx].getStridedList(length, weights, stride);
255  }
256 
257 
258  template <typename Adapter>
259  void applyPartitioningSolution(const User &in, User *&out,
260  const PartitioningSolution<Adapter> &solution) const;
261 
262  template <typename Adapter>
263  void applyPartitioningSolution(const User &in, RCP<User> &out,
264  const PartitioningSolution<Adapter> &solution) const;
265 
266 private:
267 
268  RCP<const User > ingraph_;
269  RCP<const xgraph_t > graph_;
270  RCP<const Comm<int> > comm_;
271 
272  ArrayRCP<const offset_t> offs_;
273  ArrayRCP<const gno_t> adjids_;
274 
275  int nWeightsPerVertex_;
276  ArrayRCP<StridedData<lno_t, scalar_t> > vertexWeights_;
277  ArrayRCP<bool> vertexDegreeWeight_;
278 
279  int nWeightsPerEdge_;
280  ArrayRCP<StridedData<lno_t, scalar_t> > edgeWeights_;
281 
282  int coordinateDim_;
283  ArrayRCP<StridedData<lno_t, scalar_t> > coords_;
284 
285 };
286 
288 // Definitions
290 
291 template <typename User, typename UserCoord>
293  const RCP<const User> &ingraph, int nVtxWgts, int nEdgeWgts):
294  ingraph_(ingraph), graph_(), comm_() , offs_(), adjids_(),
295  nWeightsPerVertex_(nVtxWgts), vertexWeights_(), vertexDegreeWeight_(),
296  nWeightsPerEdge_(nEdgeWgts), edgeWeights_(),
297  coordinateDim_(0), coords_()
298 {
299  typedef StridedData<lno_t,scalar_t> input_t;
300 
301  try {
302  graph_ = rcp_const_cast<const xgraph_t>(
303  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(ingraph)));
304  }
306 
307  comm_ = graph_->getComm();
308  size_t nvtx = graph_->getNodeNumRows();
309  size_t nedges = graph_->getNodeNumEntries();
310 
311  // Unfortunately we have to copy the offsets and edge Ids
312  // because edge Ids are not usually stored in vertex id order.
313 
314  size_t n = nvtx + 1;
315  offset_t *offs = new offset_t [n];
316 
317  if (!offs)
318  {
319  std::cerr << "Error: " << __FILE__ << ", " << __LINE__<< std::endl;
320  std::cerr << n << " objects" << std::endl;
321  throw std::bad_alloc();
322  }
323 
324  gno_t *adjids = NULL;
325  if (nedges)
326  {
327  adjids = new gno_t [nedges];
328 
329  if (!adjids)
330  {
331  std::cerr << "Error: " << __FILE__ << ", " << __LINE__<< std::endl;
332  std::cerr << nedges << " objects" << std::endl;
333  throw std::bad_alloc();
334  }
335  }
336 
337  offs[0] = 0;
338  for (size_t v=0; v < nvtx; v++){
339  ArrayView<const lno_t> nbors;
340  graph_->getLocalRowView(v, nbors);
341  offs[v+1] = offs[v] + nbors.size();
342  for (offset_t e=offs[v], i=0; e < offs[v+1]; e++)
343  adjids[e] = graph_->getColMap()->getGlobalElement(nbors[i++]);
344  }
345 
346  offs_ = arcp(offs, 0, n, true);
347  adjids_ = arcp(adjids, 0, nedges, true);
348 
349  if (nWeightsPerVertex_ > 0) {
350  vertexWeights_ =
351  arcp(new input_t[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
352  vertexDegreeWeight_ =
353  arcp(new bool[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
354  for (int i=0; i < nWeightsPerVertex_; i++)
355  vertexDegreeWeight_[i] = false;
356  }
357 
358  if (nWeightsPerEdge_ > 0)
359  edgeWeights_ = arcp(new input_t[nWeightsPerEdge_], 0, nWeightsPerEdge_, true);
360 }
361 
363 template <typename User, typename UserCoord>
365  const scalar_t *weightVal, int stride, int idx)
366 {
367  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
368  setVertexWeights(weightVal, stride, idx);
369  else
370  setEdgeWeights(weightVal, stride, idx);
371 }
372 
374 template <typename User, typename UserCoord>
376  const scalar_t *weightVal, int stride, int idx)
377 {
378  typedef StridedData<lno_t,scalar_t> input_t;
379 
380  if(idx<0 || idx >= nWeightsPerVertex_)
381  {
382  std::ostringstream emsg;
383  emsg << __FILE__ << ":" << __LINE__
384  << " Invalid vertex weight index " << idx << std::endl;
385  throw std::runtime_error(emsg.str());
386  }
387 
388  size_t nvtx = getLocalNumVertices();
389  ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride, false);
390  vertexWeights_[idx] = input_t(weightV, stride);
391 }
392 
394 template <typename User, typename UserCoord>
396  int idx)
397 {
398  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
399  setVertexWeightIsDegree(idx);
400  else {
401  std::ostringstream emsg;
402  emsg << __FILE__ << "," << __LINE__
403  << " error: setWeightIsNumberOfNonZeros is supported only for"
404  << " vertices" << std::endl;
405  throw std::runtime_error(emsg.str());
406  }
407 }
408 
410 template <typename User, typename UserCoord>
412  int idx)
413 {
414  if(idx<0 || idx >= nWeightsPerVertex_)
415  {
416  std::ostringstream emsg;
417  emsg << __FILE__ << ":" << __LINE__
418  << " Invalid vertex weight index " << idx << std::endl;
419  throw std::runtime_error(emsg.str());
420  }
421 
422  vertexDegreeWeight_[idx] = true;
423 }
424 
426 template <typename User, typename UserCoord>
428  const scalar_t *weightVal, int stride, int idx)
429 {
430  typedef StridedData<lno_t,scalar_t> input_t;
431 
432  if(idx<0 || idx >= nWeightsPerEdge_)
433  {
434  std::ostringstream emsg;
435  emsg << __FILE__ << ":" << __LINE__
436  << " Invalid edge weight index " << idx << std::endl;
437  throw std::runtime_error(emsg.str());
438  }
439 
440  size_t nedges = getLocalNumEdges();
441  ArrayRCP<const scalar_t> weightV(weightVal, 0, nedges*stride, false);
442  edgeWeights_[idx] = input_t(weightV, stride);
443 }
444 
446 template <typename User, typename UserCoord>
447  template<typename Adapter>
449  const User &in, User *&out,
450  const PartitioningSolution<Adapter> &solution) const
451 {
452  // Get an import list (rows to be received)
453  size_t numNewVtx;
454  ArrayRCP<gno_t> importList;
455  try{
456  numNewVtx = Zoltan2::getImportList<Adapter,
458  (solution, this, importList);
459  }
461 
462  // Move the rows, creating a new graph.
463  RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewVtx,
464  importList.getRawPtr());
465  out = outPtr.get();
466  outPtr.release();
467 }
468 
470 template <typename User, typename UserCoord>
471  template<typename Adapter>
473  const User &in, RCP<User> &out,
474  const PartitioningSolution<Adapter> &solution) const
475 {
476  // Get an import list (rows to be received)
477  size_t numNewVtx;
478  ArrayRCP<gno_t> importList;
479  try{
480  numNewVtx = Zoltan2::getImportList<Adapter,
482  (solution, this, importList);
483  }
485 
486  // Move the rows, creating a new graph.
487  out = XpetraTraits<User>::doMigration(in, numNewVtx,
488  importList.getRawPtr());
489 }
490 
491 } //namespace Zoltan2
492 
493 #endif
InputTraits< User >::scalar_t scalar_t
bool useDegreeAsVertexWeight(int idx) const
Indicate whether vertex weight with index idx should be the global degree of the vertex.
Helper functions for Partitioning Problems.
RCP< const xgraph_t > getXpetraGraph() const
Access to Xpetra-wrapped user&#39;s graph.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
InputTraits< User >::gno_t gno_t
GraphAdapter defines the interface for graph-based user data.
void setWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to weights for the primary entity type.
default_part_t part_t
The data type to represent part numbers.
default_offset_t offset_t
The data type to represent offsets.
InputTraits< User >::lno_t lno_t
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...
Provides access for Zoltan2 to Xpetra::CrsGraph data.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
Traits of Xpetra classes, including migration method.
int getNumWeightsPerEdge() const
Returns the number (0 or greater) of edge weights.
RCP< const User > getUserGraph() const
Access to user&#39;s graph.
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
size_t getLocalNumEdges() const
Returns the number of edges on this process.
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
XpetraCrsGraphAdapter(const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
A PartitioningSolution is a solution to a partitioning problem.
int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
The StridedData class manages lists of weights or coordinates.
InputTraits< User >::part_t part_t
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.
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
void setEdgeWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to edge weights.
void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
void getVertexIDsView(const gno_t *&ids) const
Sets pointers to this process&#39; graph entries.
Defines the GraphAdapter interface.
void getVertexWeightsView(const scalar_t *&weights, int &stride, int idx) const
void setVertexWeightIsDegree(int idx)
Specify an index for which the vertex weight should be the degree of the vertex.
InputTraits< User >::offset_t offset_t
void setVertexWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to vertex weights.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
default_scalar_t scalar_t
The data type for weights and coordinates.
This file defines the StridedData class.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74
size_t getLocalNumVertices() const
Returns the number of vertices on this process.