Zoltan2
Zoltan2_TpetraRowGraphAdapter.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_TPETRAROWGRAPHADAPTER_HPP_
51 #define _ZOLTAN2_TPETRAROWGRAPHADAPTER_HPP_
52 
53 #include <Zoltan2_GraphAdapter.hpp>
54 #include <Zoltan2_StridedData.hpp>
56 #include <Tpetra_RowGraph.hpp>
57 
58 namespace Zoltan2 {
59 
81 template <typename User, typename UserCoord=User>
82  class TpetraRowGraphAdapter : public GraphAdapter<User,UserCoord> {
83 
84 public:
85 
86 #ifndef DOXYGEN_SHOULD_SKIP_THIS
87  typedef typename InputTraits<User>::scalar_t scalar_t;
88  typedef typename InputTraits<User>::lno_t lno_t;
89  typedef typename InputTraits<User>::gno_t gno_t;
90  typedef typename InputTraits<User>::part_t part_t;
91  typedef typename InputTraits<User>::node_t node_t;
92  typedef User user_t;
93  typedef UserCoord userCoord_t;
94 #endif
95 
99 
109  TpetraRowGraphAdapter(const RCP<const User> &ingraph,
110  int nVtxWeights=0, int nEdgeWeights=0);
111 
124  void setWeights(const scalar_t *val, int stride, int idx);
125 
141  void setVertexWeights(const scalar_t *val, int stride, int idx);
142 
148  void setWeightIsDegree(int idx);
149 
155  void setVertexWeightIsDegree(int idx);
156 
179  void setEdgeWeights(const scalar_t *val, int stride, int idx);
180 
182  // The Adapter interface.
184 
186  // The GraphAdapter interface.
188 
189  // TODO: Assuming rows == objects;
190  // TODO: Need to add option for columns or nonzeros?
191  size_t getLocalNumVertices() const { return graph_->getNodeNumRows(); }
192 
193  void getVertexIDsView(const gno_t *&ids) const
194  {
195  ids = NULL;
196  if (getLocalNumVertices())
197  ids = graph_->getRowMap()->getNodeElementList().getRawPtr();
198  }
199 
200  size_t getLocalNumEdges() const { return graph_->getNodeNumEntries(); }
201 
202  void getEdgesView(const lno_t *&offsets, const gno_t *&adjIds) const
203  {
204  offsets = offs_.getRawPtr();
205  adjIds = (getLocalNumEdges() ? adjids_.getRawPtr() : NULL);
206  }
207 
208  int getNumWeightsPerVertex() const { return nWeightsPerVertex_;}
209 
210  void getVertexWeightsView(const scalar_t *&weights, int &stride,
211  int idx) const
212  {
213  env_->localInputAssertion(__FILE__, __LINE__, "invalid weight index",
214  idx >= 0 && idx < nWeightsPerVertex_, BASIC_ASSERTION);
215  size_t length;
216  vertexWeights_[idx].getStridedList(length, weights, stride);
217  }
218 
219  bool useDegreeAsVertexWeight(int idx) const {return vertexDegreeWeight_[idx];}
220 
221  int getNumWeightsPerEdge() const { return nWeightsPerEdge_;}
222 
223  void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
224  {
225  env_->localInputAssertion(__FILE__, __LINE__, "invalid weight index",
226  idx >= 0 && idx < nWeightsPerEdge_, BASIC_ASSERTION);
227  size_t length;
228  edgeWeights_[idx].getStridedList(length, weights, stride);
229  }
230 
231 
232  template <typename Adapter>
233  void applyPartitioningSolution(const User &in, User *&out,
234  const PartitioningSolution<Adapter> &solution) const;
235 
236  template <typename Adapter>
237  void applyPartitioningSolution(const User &in, RCP<User> &out,
238  const PartitioningSolution<Adapter> &solution) const;
239 
240 private:
241 
242  RCP<const Environment> env_; // for error messages, etc.
243 
244  RCP<const User> graph_;
245 
246  ArrayRCP<const lno_t> offs_;
247  ArrayRCP<const gno_t> adjids_;
248 
249  int nWeightsPerVertex_;
250  ArrayRCP<StridedData<lno_t, scalar_t> > vertexWeights_;
251  ArrayRCP<bool> vertexDegreeWeight_;
252 
253  int nWeightsPerEdge_;
254  ArrayRCP<StridedData<lno_t, scalar_t> > edgeWeights_;
255 
256  int coordinateDim_;
257  ArrayRCP<StridedData<lno_t, scalar_t> > coords_;
258 
259  RCP<User> doMigration(const User &from, size_t numLocalRows,
260  const gno_t *myNewRows) const;
261 };
262 
264 // Definitions
266 
267 template <typename User, typename UserCoord>
269  const RCP<const User> &ingraph, int nVtxWgts, int nEdgeWgts):
270  env_(rcp(new Environment)),
271  graph_(ingraph), offs_(),
272  adjids_(),
273  nWeightsPerVertex_(nVtxWgts), vertexWeights_(), vertexDegreeWeight_(),
274  nWeightsPerEdge_(nEdgeWgts), edgeWeights_(),
275  coordinateDim_(0), coords_()
276 {
277  typedef StridedData<lno_t,scalar_t> input_t;
278 
279  size_t nvtx = graph_->getNodeNumRows();
280  size_t nedges = graph_->getNodeNumEntries();
281  size_t maxnumentries =
282  graph_->getNodeMaxNumRowEntries(); // Diff from CrsMatrix
283 
284  // Unfortunately we have to copy the offsets and edge Ids
285  // because edge Ids are not usually stored in vertex id order.
286 
287  size_t n = nvtx + 1;
288  lno_t *offs = new lno_t [n];
289  env_->localMemoryAssertion(__FILE__, __LINE__, n, offs);
290 
291  gno_t *adjids = NULL;
292  if (nedges){
293  adjids = new gno_t [nedges];
294  env_->localMemoryAssertion(__FILE__, __LINE__, nedges, adjids);
295  }
296 
297  ArrayRCP<lno_t> nbors(maxnumentries); // Diff from CrsGraph
298 
299  offs[0] = 0;
300  for (size_t v=0; v < nvtx; v++){
301  graph_->getLocalRowCopy(v, nbors(), nedges); // Diff from CrsGraph
302  offs[v+1] = offs[v] + nedges;
303  for (lno_t e=offs[v], i=0; e < offs[v+1]; e++) {
304  adjids[e] = graph_->getColMap()->getGlobalElement(nbors[i++]);
305  }
306  }
307 
308  offs_ = arcp(offs, 0, n, true);
309  adjids_ = arcp(adjids, 0, nedges, true);
310 
311  if (nWeightsPerVertex_ > 0) {
312  vertexWeights_ =
313  arcp(new input_t[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
314  vertexDegreeWeight_ =
315  arcp(new bool[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
316  for (int i=0; i < nWeightsPerVertex_; i++)
317  vertexDegreeWeight_[i] = false;
318  }
319 
320  if (nWeightsPerEdge_ > 0)
321  edgeWeights_ = arcp(new input_t[nWeightsPerEdge_], 0, nWeightsPerEdge_, true);
322 }
323 
325 template <typename User, typename UserCoord>
327  const scalar_t *weightVal, int stride, int idx)
328 {
329  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
330  setVertexWeights(weightVal, stride, idx);
331  else
332  setEdgeWeights(weightVal, stride, idx);
333 }
334 
336 template <typename User, typename UserCoord>
338  const scalar_t *weightVal, int stride, int idx)
339 {
340  typedef StridedData<lno_t,scalar_t> input_t;
341  env_->localInputAssertion(__FILE__, __LINE__, "invalid vertex weight index",
342  idx >= 0 && idx < nWeightsPerVertex_, BASIC_ASSERTION);
343  size_t nvtx = getLocalNumVertices();
344  ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride, false);
345  vertexWeights_[idx] = input_t(weightV, stride);
346 }
347 
349 template <typename User, typename UserCoord>
351  int idx)
352 {
353  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
355  else {
356  std::ostringstream emsg;
357  emsg << __FILE__ << "," << __LINE__
358  << " error: setWeightIsNumberOfNonZeros is supported only for"
359  << " vertices" << std::endl;
360  throw std::runtime_error(emsg.str());
361  }
362 }
363 
365 template <typename User, typename UserCoord>
367  int idx)
368 {
369  env_->localInputAssertion(__FILE__, __LINE__, "invalid vertex weight index",
370  idx >= 0 && idx < nWeightsPerVertex_, BASIC_ASSERTION);
371 
372  vertexDegreeWeight_[idx] = true;
373 }
374 
376 template <typename User, typename UserCoord>
378  const scalar_t *weightVal, int stride, int idx)
379 {
380  typedef StridedData<lno_t,scalar_t> input_t;
381  env_->localInputAssertion(__FILE__, __LINE__, "invalid edge weight index",
382  idx >= 0 && idx < nWeightsPerEdge_, BASIC_ASSERTION);
383  size_t nedges = getLocalNumEdges();
384  ArrayRCP<const scalar_t> weightV(weightVal, 0, nedges*stride, false);
385  edgeWeights_[idx] = input_t(weightV, stride);
386 }
387 
389 template <typename User, typename UserCoord>
390  template<typename Adapter>
392  const User &in, User *&out,
393  const PartitioningSolution<Adapter> &solution) const
394 {
395  // Get an import list (rows to be received)
396  size_t numNewVtx;
397  ArrayRCP<gno_t> importList;
398  try{
399  numNewVtx = Zoltan2::getImportList<Adapter,
401  (solution, this, importList);
402  }
404 
405  // Move the rows, creating a new graph.
406  RCP<User> outPtr = doMigration(in, numNewVtx, importList.getRawPtr());
407  out = outPtr.get();
408  outPtr.release();
409 }
410 
412 template <typename User, typename UserCoord>
413  template<typename Adapter>
415  const User &in, RCP<User> &out,
416  const PartitioningSolution<Adapter> &solution) const
417 {
418  // Get an import list (rows to be received)
419  size_t numNewVtx;
420  ArrayRCP<gno_t> importList;
421  try{
422  numNewVtx = Zoltan2::getImportList<Adapter,
424  (solution, this, importList);
425  }
427 
428  // Move the rows, creating a new graph.
429  out = doMigration(in, numNewVtx, importList.getRawPtr());
430 }
431 
432 
434 template < typename User, typename UserCoord>
436  const User &from,
437  size_t numLocalRows,
438  const gno_t *myNewRows
439 ) const
440 {
441  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
442  typedef Tpetra::CrsGraph<lno_t, gno_t, node_t> tcrsgraph_t;
443  typedef Tpetra::RowGraph<lno_t, gno_t, node_t> trowgraph_t;
444 
445  // We cannot create a Tpetra::RowGraph, unless the underlying type is
446  // something we know (like Tpetra::CrsGraph).
447  // If the underlying type is something different, the user probably doesn't
448  // want a Tpetra::CrsGraph back, so we throw an error.
449 
450  // Try to cast "from" graph to a TPetra::CrsGraph
451  // If that fails we throw an error.
452  // We could cast as a ref which will throw std::bad_cast but with ptr
453  // approach it might be clearer what's going on here
454  const tcrsgraph_t *pCrsGraphSrc = dynamic_cast<const tcrsgraph_t *>(&from);
455 
456  if(!pCrsGraphSrc) {
457  throw std::logic_error("TpetraRowGraphAdapter cannot migrate data for "
458  "your RowGraph; it can migrate data only for "
459  "Tpetra::CrsGraph. "
460  "You can inherit from TpetraRowGraphAdapter and "
461  "implement migration for your RowGraph.");
462  }
463 
464  lno_t base = 0;
465 
466  // source map
467  const RCP<const map_t> &smap = from.getRowMap();
468  int oldNumElts = smap->getNodeNumElements();
469  gno_t numGlobalRows = smap->getGlobalNumElements();
470 
471  // target map
472  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
473  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
474  RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
475 
476  // importer
477  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
478 
479  // number of entries in my new rows
480  typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
481  vector_t numOld(smap);
482  vector_t numNew(tmap);
483  for (int lid=0; lid < oldNumElts; lid++){
484  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
485  from.getNumEntriesInLocalRow(lid));
486  }
487  numNew.doImport(numOld, importer, Tpetra::INSERT);
488 
489  size_t numElts = tmap->getNodeNumElements();
490  ArrayRCP<const gno_t> nnz;
491  if (numElts > 0)
492  nnz = numNew.getData(0); // hangs if vector len == 0
493 
494  ArrayRCP<const size_t> nnz_size_t;
495 
496  if (numElts && sizeof(gno_t) != sizeof(size_t)){
497  size_t *vals = new size_t [numElts];
498  nnz_size_t = arcp(vals, 0, numElts, true);
499  for (size_t i=0; i < numElts; i++){
500  vals[i] = static_cast<size_t>(nnz[i]);
501  }
502  }
503  else{
504  nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
505  }
506 
507  // target graph
508  RCP<tcrsgraph_t> G = rcp(new tcrsgraph_t(tmap, nnz_size_t,
509  Tpetra::StaticProfile));
510 
511  G->doImport(*pCrsGraphSrc, importer, Tpetra::INSERT);
512  G->fillComplete();
513  return Teuchos::rcp_dynamic_cast<User>(G);
514 }
515 
516 } //namespace Zoltan2
517 
518 #endif
InputTraits< User >::scalar_t scalar_t
Helper functions for Partitioning Problems.
fast typical checks for valid arguments
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
size_t getLocalNumEdges() const
Returns the number of edges on this process.
InputTraits< User >::gno_t gno_t
bool useDegreeAsVertexWeight(int idx) const
Indicate whether vertex weight with index idx should be the global degree of the vertex.
GraphAdapter defines the interface for graph-based user data.
void setVertexWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to vertex weights.
default_part_t part_t
The data type to represent part numbers.
InputTraits< User >::lno_t lno_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights
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...
void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
int getNumWeightsPerEdge() const
Returns the number (0 or greater) of edge weights.
dictionary vals
Definition: xml2dox.py:186
A PartitioningSolution is a solution to a partitioning problem.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
The StridedData class manages lists of weights or coordinates.
InputTraits< User >::part_t part_t
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
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.
size_t getLocalNumVertices() const
Returns the number of vertices 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 setVertexWeightIsDegree(int idx)
Specify an index for which the vertex weight should be the degree of the vertex.
TpetraRowGraphAdapter(const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
void setEdgeWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to edge weights.
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 setWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to weights for the primary entity type.
default_scalar_t scalar_t
The data type for weights and coordinates.
void getEdgesView(const lno_t *&offsets, const gno_t *&adjIds) const
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
Provides access for Zoltan2 to Tpetra::RowGraph data.
This file defines the StridedData class.