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>::offset_t offset_t;
89  typedef typename InputTraits<User>::lno_t lno_t;
90  typedef typename InputTraits<User>::gno_t gno_t;
91  typedef typename InputTraits<User>::part_t part_t;
92  typedef typename InputTraits<User>::node_t node_t;
93  typedef User user_t;
94  typedef UserCoord userCoord_t;
95 #endif
96 
100 
110  TpetraRowGraphAdapter(const RCP<const User> &ingraph,
111  int nVtxWeights=0, int nEdgeWeights=0);
112 
125  void setWeights(const scalar_t *val, int stride, int idx);
126 
142  void setVertexWeights(const scalar_t *val, int stride, int idx);
143 
149  void setWeightIsDegree(int idx);
150 
156  void setVertexWeightIsDegree(int idx);
157 
180  void setEdgeWeights(const scalar_t *val, int stride, int idx);
181 
183  // The Adapter interface.
185 
187  // The GraphAdapter interface.
189 
190  // TODO: Assuming rows == objects;
191  // TODO: Need to add option for columns or nonzeros?
192  size_t getLocalNumVertices() const { return graph_->getNodeNumRows(); }
193 
194  void getVertexIDsView(const gno_t *&ids) const
195  {
196  ids = NULL;
197  if (getLocalNumVertices())
198  ids = graph_->getRowMap()->getNodeElementList().getRawPtr();
199  }
200 
201  size_t getLocalNumEdges() const { return graph_->getNodeNumEntries(); }
202 
203  void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
204  {
205  offsets = offs_.getRawPtr();
206  adjIds = (getLocalNumEdges() ? adjids_.getRawPtr() : NULL);
207  }
208 
209  int getNumWeightsPerVertex() const { return nWeightsPerVertex_;}
210 
211  void getVertexWeightsView(const scalar_t *&weights, int &stride,
212  int idx) const
213  {
214  if(idx<0 || idx >= nWeightsPerVertex_)
215  {
216  std::ostringstream emsg;
217  emsg << __FILE__ << ":" << __LINE__
218  << " Invalid vertex weight index " << idx << std::endl;
219  throw std::runtime_error(emsg.str());
220  }
221 
222  size_t length;
223  vertexWeights_[idx].getStridedList(length, weights, stride);
224  }
225 
226  bool useDegreeAsVertexWeight(int idx) const {return vertexDegreeWeight_[idx];}
227 
228  int getNumWeightsPerEdge() const { return nWeightsPerEdge_;}
229 
230  void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
231  {
232  if(idx<0 || idx >= nWeightsPerEdge_)
233  {
234  std::ostringstream emsg;
235  emsg << __FILE__ << ":" << __LINE__
236  << " Invalid edge weight index " << idx << std::endl;
237  throw std::runtime_error(emsg.str());
238  }
239 
240  size_t length;
241  edgeWeights_[idx].getStridedList(length, weights, stride);
242  }
243 
244 
245  template <typename Adapter>
246  void applyPartitioningSolution(const User &in, User *&out,
247  const PartitioningSolution<Adapter> &solution) const;
248 
249  template <typename Adapter>
250  void applyPartitioningSolution(const User &in, RCP<User> &out,
251  const PartitioningSolution<Adapter> &solution) const;
252 
253 private:
254 
255  RCP<const User> graph_;
256 
257  ArrayRCP<const offset_t> offs_;
258  ArrayRCP<const gno_t> adjids_;
259 
260  int nWeightsPerVertex_;
261  ArrayRCP<StridedData<lno_t, scalar_t> > vertexWeights_;
262  ArrayRCP<bool> vertexDegreeWeight_;
263 
264  int nWeightsPerEdge_;
265  ArrayRCP<StridedData<lno_t, scalar_t> > edgeWeights_;
266 
267  int coordinateDim_;
268  ArrayRCP<StridedData<lno_t, scalar_t> > coords_;
269 
270  RCP<User> doMigration(const User &from, size_t numLocalRows,
271  const gno_t *myNewRows) const;
272 };
273 
275 // Definitions
277 
278 template <typename User, typename UserCoord>
280  const RCP<const User> &ingraph, int nVtxWgts, int nEdgeWgts):
281  graph_(ingraph), offs_(),
282  adjids_(),
283  nWeightsPerVertex_(nVtxWgts), vertexWeights_(), vertexDegreeWeight_(),
284  nWeightsPerEdge_(nEdgeWgts), edgeWeights_(),
285  coordinateDim_(0), coords_()
286 {
287  typedef StridedData<lno_t,scalar_t> input_t;
288 
289  size_t nvtx = graph_->getNodeNumRows();
290  size_t nedges = graph_->getNodeNumEntries();
291  size_t maxnumentries =
292  graph_->getNodeMaxNumRowEntries(); // Diff from CrsMatrix
293 
294  // Unfortunately we have to copy the offsets and edge Ids
295  // because edge Ids are not usually stored in vertex id order.
296 
297  size_t n = nvtx + 1;
298  offset_t *offs = new offset_t [n];
299 
300  if (!offs)
301  {
302  std::cerr << "Error: " << __FILE__ << ", " << __LINE__<< std::endl;
303  std::cerr << n << " objects" << std::endl;
304  throw std::bad_alloc();
305  }
306 
307  gno_t *adjids = NULL;
308  if (nedges)
309  {
310  adjids = new gno_t [nedges];
311 
312  if (!adjids)
313  {
314  std::cerr << "Error: " << __FILE__ << ", " << __LINE__<< std::endl;
315  std::cerr << nedges << " objects" << std::endl;
316  throw std::bad_alloc();
317  }
318  }
319 
320  typename User::nonconst_local_inds_host_view_type nbors("nbors", maxnumentries);
321 
322  offs[0] = 0;
323  for (size_t v=0; v < nvtx; v++){
324  graph_->getLocalRowCopy(v, nbors, nedges); // Diff from CrsGraph
325  offs[v+1] = offs[v] + nedges;
326  for (offset_t e=offs[v], i=0; e < offs[v+1]; e++) {
327  adjids[e] = graph_->getColMap()->getGlobalElement(nbors[i++]);
328  }
329  }
330 
331  offs_ = arcp(offs, 0, n, true);
332  adjids_ = arcp(adjids, 0, nedges, true);
333 
334  if (nWeightsPerVertex_ > 0) {
335  vertexWeights_ =
336  arcp(new input_t[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
337  vertexDegreeWeight_ =
338  arcp(new bool[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
339  for (int i=0; i < nWeightsPerVertex_; i++)
340  vertexDegreeWeight_[i] = false;
341  }
342 
343  if (nWeightsPerEdge_ > 0)
344  edgeWeights_ = arcp(new input_t[nWeightsPerEdge_], 0, nWeightsPerEdge_, true);
345 }
346 
348 template <typename User, typename UserCoord>
350  const scalar_t *weightVal, int stride, int idx)
351 {
352  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
353  setVertexWeights(weightVal, stride, idx);
354  else
355  setEdgeWeights(weightVal, stride, idx);
356 }
357 
359 template <typename User, typename UserCoord>
361  const scalar_t *weightVal, int stride, int idx)
362 {
363  typedef StridedData<lno_t,scalar_t> input_t;
364  if(idx<0 || idx >= nWeightsPerVertex_)
365  {
366  std::ostringstream emsg;
367  emsg << __FILE__ << ":" << __LINE__
368  << " Invalid vertex weight index " << idx << std::endl;
369  throw std::runtime_error(emsg.str());
370  }
371 
372  size_t nvtx = getLocalNumVertices();
373  ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride, false);
374  vertexWeights_[idx] = input_t(weightV, stride);
375 }
376 
378 template <typename User, typename UserCoord>
380  int idx)
381 {
382  if (this->getPrimaryEntityType() == GRAPH_VERTEX)
383  setVertexWeightIsDegree(idx);
384  else {
385  std::ostringstream emsg;
386  emsg << __FILE__ << "," << __LINE__
387  << " error: setWeightIsNumberOfNonZeros is supported only for"
388  << " vertices" << std::endl;
389  throw std::runtime_error(emsg.str());
390  }
391 }
392 
394 template <typename User, typename UserCoord>
396  int idx)
397 {
398  if(idx<0 || idx >= nWeightsPerVertex_)
399  {
400  std::ostringstream emsg;
401  emsg << __FILE__ << ":" << __LINE__
402  << " Invalid vertex weight index " << idx << std::endl;
403  throw std::runtime_error(emsg.str());
404  }
405 
406  vertexDegreeWeight_[idx] = true;
407 }
408 
410 template <typename User, typename UserCoord>
412  const scalar_t *weightVal, int stride, int idx)
413 {
414  typedef StridedData<lno_t,scalar_t> input_t;
415 
416  if(idx<0 || idx >= nWeightsPerEdge_)
417  {
418  std::ostringstream emsg;
419  emsg << __FILE__ << ":" << __LINE__
420  << " Invalid edge weight index " << idx << std::endl;
421  throw std::runtime_error(emsg.str());
422  }
423 
424  size_t nedges = getLocalNumEdges();
425  ArrayRCP<const scalar_t> weightV(weightVal, 0, nedges*stride, false);
426  edgeWeights_[idx] = input_t(weightV, stride);
427 }
428 
430 template <typename User, typename UserCoord>
431  template<typename Adapter>
433  const User &in, User *&out,
434  const PartitioningSolution<Adapter> &solution) const
435 {
436  // Get an import list (rows to be received)
437  size_t numNewVtx;
438  ArrayRCP<gno_t> importList;
439  try{
440  numNewVtx = Zoltan2::getImportList<Adapter,
442  (solution, this, importList);
443  }
445 
446  // Move the rows, creating a new graph.
447  RCP<User> outPtr = doMigration(in, numNewVtx, importList.getRawPtr());
448  out = outPtr.get();
449  outPtr.release();
450 }
451 
453 template <typename User, typename UserCoord>
454  template<typename Adapter>
456  const User &in, RCP<User> &out,
457  const PartitioningSolution<Adapter> &solution) const
458 {
459  // Get an import list (rows to be received)
460  size_t numNewVtx;
461  ArrayRCP<gno_t> importList;
462  try{
463  numNewVtx = Zoltan2::getImportList<Adapter,
465  (solution, this, importList);
466  }
468 
469  // Move the rows, creating a new graph.
470  out = doMigration(in, numNewVtx, importList.getRawPtr());
471 }
472 
473 
475 template < typename User, typename UserCoord>
477  const User &from,
478  size_t numLocalRows,
479  const gno_t *myNewRows
480 ) const
481 {
482  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
483  typedef Tpetra::CrsGraph<lno_t, gno_t, node_t> tcrsgraph_t;
484 
485  // We cannot create a Tpetra::RowGraph, unless the underlying type is
486  // something we know (like Tpetra::CrsGraph).
487  // If the underlying type is something different, the user probably doesn't
488  // want a Tpetra::CrsGraph back, so we throw an error.
489 
490  // Try to cast "from" graph to a TPetra::CrsGraph
491  // If that fails we throw an error.
492  // We could cast as a ref which will throw std::bad_cast but with ptr
493  // approach it might be clearer what's going on here
494  const tcrsgraph_t *pCrsGraphSrc = dynamic_cast<const tcrsgraph_t *>(&from);
495 
496  if(!pCrsGraphSrc) {
497  throw std::logic_error("TpetraRowGraphAdapter cannot migrate data for "
498  "your RowGraph; it can migrate data only for "
499  "Tpetra::CrsGraph. "
500  "You can inherit from TpetraRowGraphAdapter and "
501  "implement migration for your RowGraph.");
502  }
503 
504  // source map
505  const RCP<const map_t> &smap = from.getRowMap();
506  int oldNumElts = smap->getNodeNumElements();
507  gno_t numGlobalRows = smap->getGlobalNumElements();
508  gno_t base = smap->getMinAllGlobalIndex();
509 
510  // target map
511  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
512  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
513  RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
514 
515  // importer
516  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
517 
518  // number of entries in my new rows
519  typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
520  vector_t numOld(smap);
521  vector_t numNew(tmap);
522  for (int lid=0; lid < oldNumElts; lid++){
523  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
524  from.getNumEntriesInLocalRow(lid));
525  }
526  numNew.doImport(numOld, importer, Tpetra::INSERT);
527 
528  size_t numElts = tmap->getNodeNumElements();
529  ArrayRCP<const gno_t> nnz;
530  if (numElts > 0)
531  nnz = numNew.getData(0); // hangs if vector len == 0
532 
533  ArrayRCP<const size_t> nnz_size_t;
534 
535  if (numElts && sizeof(gno_t) != sizeof(size_t)){
536  size_t *vals = new size_t [numElts];
537  nnz_size_t = arcp(vals, 0, numElts, true);
538  for (size_t i=0; i < numElts; i++){
539  vals[i] = static_cast<size_t>(nnz[i]);
540  }
541  }
542  else{
543  nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
544  }
545 
546  // target graph
547  RCP<tcrsgraph_t> G =
548  rcp(new tcrsgraph_t(tmap, nnz_size_t(), Tpetra::StaticProfile));
549 
550  G->doImport(*pCrsGraphSrc, importer, Tpetra::INSERT);
551  G->fillComplete();
552  return Teuchos::rcp_dynamic_cast<User>(G);
553 }
554 
555 } //namespace Zoltan2
556 
557 #endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the GraphAdapter interface.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
InputTraits< User >::node_t node_t
InputTraits< User >::offset_t offset_t
InputTraits< User >::part_t part_t
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
InputTraits< User >::gno_t gno_t
GraphAdapter defines the interface for graph-based user data.
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
Provides access for Zoltan2 to Tpetra::RowGraph data.
size_t getLocalNumEdges() const
Returns the number of edges on this process.
void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide a pointer to the edge weights, if any.
size_t getLocalNumVertices() const
Returns the number of vertices on this process.
void getVertexIDsView(const gno_t *&ids) const
bool useDegreeAsVertexWeight(int idx) const
Indicate whether vertex weight with index idx should be the global degree of the vertex.
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
void setEdgeWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to edge weights.
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
void setWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to weights for the primary entity type.
void setVertexWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to vertex weights.
TpetraRowGraphAdapter(const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
void setVertexWeightIsDegree(int idx)
Specify an index for which the vertex weight should be the degree of the vertex.
void getVertexWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide a pointer to the vertex weights, if any.
int getNumWeightsPerEdge() const
Returns the number (0 or greater) of edge weights.
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Tpetra::Map map_t
Definition: mapRemotes.cpp:16
Created by mbenlioglu on Aug 31, 2020.
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...
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_part_t part_t
The data type to represent part numbers.
default_scalar_t scalar_t
The data type for weights and coordinates.