Zoltan2
Zoltan2_HyperGraphModel.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_HYPERGRAPHMODEL_HPP_
51 #define _ZOLTAN2_HYPERGRAPHMODEL_HPP_
52 
53 #include <Zoltan2_Model.hpp>
54 #include <Zoltan2_ModelHelpers.hpp>
55 #include <Zoltan2_InputTraits.hpp>
57 #include <Zoltan2_GraphAdapter.hpp>
60 #include <Zoltan2_StridedData.hpp>
61 #include <Zoltan2_MeshAdapter.hpp>
62 
63 #include <vector>
64 #include <unordered_map>
65 #include <queue>
66 #include <Teuchos_Hashtable.hpp>
67 
68 namespace Zoltan2 {
69 
77 };
78 
80 
89 template <typename Adapter>
90 class HyperGraphModel : public Model<Adapter>
91 {
92 public:
93 
94 #ifndef DOXYGEN_SHOULD_SKIP_THIS
95  typedef typename Adapter::scalar_t scalar_t;
96  typedef typename Adapter::offset_t offset_t;
97  typedef typename Adapter::gno_t gno_t;
98  typedef typename Adapter::lno_t lno_t;
99  typedef typename Adapter::node_t node_t;
100  typedef typename Adapter::user_t user_t;
101  typedef typename Adapter::userCoord_t userCoord_t;
102  typedef Tpetra::Map<lno_t, gno_t> map_t;
103  typedef StridedData<lno_t, scalar_t> input_t;
104 #endif
105 
108 
121  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
122  modelFlag_t &modelFlags, CentricView view)
123  {
124  throw std::runtime_error("Building HyperGraphModel from MatrixAdapter not implemented yet");
125  }
126 
128  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
129  modelFlag_t &modelFlags, CentricView view)
130  {
131  throw std::runtime_error("Building HyperGraphModel from GraphAdapter not implemented yet");
132  }
133 
134  HyperGraphModel(const RCP<const MeshAdapter<user_t> > &ia,
135  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
136  modelFlag_t &modelflags, CentricView view);
137 
139  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
140  modelFlag_t &flags, CentricView view)
141  {
142  throw std::runtime_error("cannot build HyperGraphModel from VectorAdapter");
143  }
144 
146  const RCP<const Environment> &env, const RCP<const Comm<int> > &comm,
147  modelFlag_t &flags, CentricView view)
148  {
149  throw std::runtime_error("cannot build HyperGraphModel from IdentifierAdapter");
150  }
151 
152 
155  CentricView getCentricView() const {return view_;}
156 
159  bool areVertexIDsUnique() const {return unique;}
160 
163  size_t getLocalNumVertices() const { return numLocalVertices_; }
164 
167  size_t getLocalNumOwnedVertices() const { return numOwnedVertices_; }
168 
171  size_t getGlobalNumVertices() const { return numGlobalVertices_; }
172 
177  size_t getLocalNumHyperEdges() const { return numLocalEdges_; }
178 
181  size_t getGlobalNumHyperEdges() const { return numGlobalEdges_; }
182 
185  size_t getLocalNumPins() const {return numLocalPins_; }
186 
189  int getNumWeightsPerVertex() const { return numWeightsPerVertex_; }
190 
193  int getNumWeightsPerHyperEdge() const { return nWeightsPerEdge_; }
194 
197  int getNumWeightesPerPin() const {return nWeightsPerPin_;}
198 
201  int getCoordinateDim() const { return vCoordDim_; }
202 
211  ArrayView<const gno_t> &Ids,
212  ArrayView<input_t> &wgts) const
213  {
214  size_t nv = gids_.size();
215  Ids = gids_(0, nv);
216  wgts = vWeights_.view(0, numWeightsPerVertex_);
217  return nv;
218  }
219 
225  size_t getVertexCoords(ArrayView<input_t> &xyz) const
226  {
227  size_t nv = gids_.size();
228  xyz = vCoords_.view(0, vCoordDim_);
229  return nv;
230  }
231 
238  size_t getOwnedList(ArrayView<bool> &isOwner) const
239  {
240  size_t nv = isOwner_.size();
241  isOwner = isOwner_(0, nv);
242  return nv;
243  }
244 
252  void getVertexMaps(Teuchos::RCP<const map_t>& copiesMap,
253  Teuchos::RCP<const map_t>& onetooneMap) const
254  {
255  copiesMap = mapWithCopies;
256  onetooneMap = oneToOneMap;
257  }
258 
266  size_t getEdgeList(
267  ArrayView<const gno_t> &Ids,
268  ArrayView<input_t> &wgts) const
269  {
270  size_t nv = edgeGids_.size();
271  Ids = edgeGids_(0, nv);
272  wgts = eWeights_.view(0, nWeightsPerEdge_);
273  return nv;
274  }
275 
288  size_t getPinList( ArrayView<const gno_t> &pinIds,
289  ArrayView<const offset_t> &offsets,
290  ArrayView<input_t> &wgts) const
291  {
292  pinIds = pinGids_(0, numLocalPins_);
293  offsets = offsets_.view(0, offsets_.size());
294  wgts = pWeights_.view(0, nWeightsPerPin_);
295  return pinGids_.size();
296  }
297 
298 
300  // The Model interface.
302 
303  size_t getLocalNumObjects() const { return numLocalVertices_; }
304 
305  size_t getGlobalNumObjects() const { return numGlobalVertices_; }
306 
307 private:
308 
309  struct GhostCell {
310  lno_t lid; //Assumes lno_t is signed (-1 corresponds to not on this process)
311  gno_t gid;
312  unsigned int dist;
313  GhostCell(lno_t l,gno_t g, unsigned int d) {lid=l;gid=g;dist=d;}
314  bool operator<(const struct GhostCell& other) const {return dist>other.dist;}
315  };
316  template <typename AdapterWithCoords>
317  void shared_GetVertexCoords(const AdapterWithCoords *ia);
318 
319 
320  const RCP<const Environment > env_;
321  const RCP<const Comm<int> > comm_;
322 
323  CentricView view_;
324  bool unique;
325  ArrayRCP<const gno_t> gids_; // vertices of input graph
326 
327  ArrayRCP<bool> isOwner_;
328 
329  int numWeightsPerVertex_;
330  ArrayRCP<input_t> vWeights_;
331 
332  int vCoordDim_;
333  ArrayRCP<input_t> vCoords_;
334 
335  ArrayRCP<const gno_t> edgeGids_;
336 
337  int nWeightsPerEdge_;
338  ArrayRCP<input_t> eWeights_;
339 
340  ArrayRCP<const gno_t> pinGids_;
341  ArrayRCP<const offset_t> offsets_;
342 
343  int nWeightsPerPin_;
344  ArrayRCP<input_t> pWeights_;
345 
346  // For convenience
347 
348  size_t numLocalVertices_;
349  size_t numOwnedVertices_;
350  size_t numGlobalVertices_;
351  size_t numLocalEdges_;
352  size_t numGlobalEdges_;
353  size_t numLocalPins_;
354 
355  // For unique mapping
356  Teuchos::RCP<const map_t> mapWithCopies;
357  Teuchos::RCP<const map_t> oneToOneMap;
358 
359  // For debugging
360  void print();
361 
362 };
363 
364 
366 
368 //TODO get the weights hyperedges
369 //GFD Do we need weights for pins too?
370 template <typename Adapter>
372  const RCP<const MeshAdapter<user_t> > &ia,
373  const RCP<const Environment> &env,
374  const RCP<const Comm<int> > &comm,
375  modelFlag_t &modelFlags,
376  CentricView view):
377  env_(env),
378  comm_(comm),
379  view_(view),
380  gids_(),
381  isOwner_(),
382  numWeightsPerVertex_(0),
383  vWeights_(),
384  vCoordDim_(0),
385  vCoords_(),
386  edgeGids_(),
387  nWeightsPerEdge_(0),
388  eWeights_(),
389  pinGids_(),
390  offsets_(),
391  nWeightsPerPin_(0),
392  pWeights_(),
393  numLocalVertices_(0),
394  numGlobalVertices_(0),
395  numLocalEdges_(0),
396  numGlobalEdges_(0),
397  numLocalPins_(0)
398 {
399  env_->timerStart(MACRO_TIMERS, "HyperGraphModel constructed from MeshAdapter");
400  //Model Type is either traditional or ghosting
401  // Traditional:
402  // vertices == ia->getPrimaryEntityType()
403  // hyperedges == ia->getAdjacencyEntityType()
404  // pins == first adjacency between primary and adjacency types
405  // Ghosting:
406  // vertices == ia->getPrimaryEntityType()
407  // hyperedges == ia->getPrimaryEntityType()
408  // pins == k layers of second adjacency from primary through second adjacency types
409  std::string model_type("traditional");
410  const Teuchos::ParameterList &pl = env->getParameters();
411  const Teuchos::ParameterEntry *pe2 = pl.getEntryPtr("hypergraph_model_type");
412  if (pe2){
413  model_type = pe2->getValue<std::string>(&model_type);
414  }
415 
416  // Get the hypergraph types from adapter
417  Zoltan2::MeshEntityType primaryEType = ia->getPrimaryEntityType();
418  Zoltan2::MeshEntityType adjacencyEType = ia->getAdjacencyEntityType();
419 
420  // Get the IDs of the primary entity type; these are hypergraph vertices
421  gno_t const *vtxIds=NULL;
422  try {
423  numLocalVertices_ = ia->getLocalNumOf(primaryEType);
424  ia->getIDsViewOf(primaryEType, vtxIds);
425  size_t maxId = *(std::max_element(vtxIds,vtxIds+numLocalVertices_));
426  reduceAll(*comm_,Teuchos::REDUCE_MAX,1,&maxId,&numGlobalVertices_);
427  // TODO: KDD 1/17 The above computation of numGlobalVertices_ is
428  // TODO: correct only when the vertices are consecutively numbered
429  // TODO: starting at ID 1. Github #1024
430  }
432 
433  gids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_, false);
434 
435  //A mapping from gids to lids for efficiency
436  std::unordered_map<gno_t,lno_t> lid_mapping;
437  for (size_t i=0;i<numLocalVertices_;i++)
438  lid_mapping[gids_[i]]=i;
439 
440  // Define owners for each hypergraph vertex using Tpetra
441  // one to one map. This defines each hypergraph vertex to
442  // one process in the case that the adapter has copied
443  // primary entity types
444  //If the mesh adapter knows the entities are unique we can optimize out the ownership
445  unique = ia->areEntityIDsUnique(ia->getPrimaryEntityType());
446  numOwnedVertices_=numLocalVertices_;
447  isOwner_ = ArrayRCP<bool>(numLocalVertices_,true);
448  if (!unique) {
449 
450  Tpetra::global_size_t numGlobalCoords =
451  Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
452  mapWithCopies = rcp(new map_t(numGlobalCoords, gids_(), 0, comm));
453  // TODO KDD 1/17 It would be better to use minimum GID rather than
454  // TODO zero in the above Tpetra::Map constructor. Github #1024
455  oneToOneMap = Tpetra::createOneToOne<lno_t, gno_t>(mapWithCopies);
456 
457  numOwnedVertices_=oneToOneMap->getNodeNumElements();
458  for (size_t i=0;i<numLocalVertices_;i++) {
459  isOwner_[i] = oneToOneMap->isNodeGlobalElement(gids_[i]);
460  }
461  }
462 
463 
464  if (model_type=="traditional") {
465  // Traditional: Get the IDs of the adjacency entity type;
466  // these are hypergraph hyperedges
467 
468  gno_t const *edgeIds=NULL;
469  try {
470  numLocalEdges_ = ia->getLocalNumOf(adjacencyEType);
471  ia->getIDsViewOf(adjacencyEType, edgeIds);
472  size_t maxId = *(std::max_element(edgeIds,edgeIds+numLocalEdges_));
473  reduceAll(*comm_,Teuchos::REDUCE_MAX,1,&maxId,&numGlobalEdges_);
474  }
476 
477  edgeGids_ = arcp<const gno_t>(edgeIds, 0, numLocalEdges_, false);
478  }
479  else if (model_type=="ghosting") {
480  // Ghosting: Use the vertices as the hyperedges as well
481  numLocalEdges_ = numLocalVertices_;
482  edgeGids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_, false);
483  numGlobalEdges_ = numGlobalVertices_;
484  }
485 
486  //Define the entity types to use for the pins based on the centric view
487  Zoltan2::MeshEntityType primaryPinType = primaryEType;
488  Zoltan2::MeshEntityType adjacencyPinType = adjacencyEType;
489  size_t numPrimaryPins = numLocalVertices_;
490  if (view_==HYPEREDGE_CENTRIC) {
491  primaryPinType = adjacencyEType;
492  adjacencyPinType = primaryEType;
493  numPrimaryPins = numLocalEdges_;
494  }
495  if (model_type=="traditional") {
496  //Get the pins from using the traditional method of first adjacency
497  gno_t const *nborIds=NULL;
498  offset_t const *offsets=NULL;
499 
500  try {
501  ia->getAdjsView(primaryPinType,adjacencyPinType,offsets,nborIds);
502  }
504 
505  numLocalPins_ = offsets[numPrimaryPins];
506 
507  pinGids_ = arcp<const gno_t>(nborIds, 0, numLocalPins_, false);
508  offsets_ = arcp<const offset_t>(offsets, 0, numPrimaryPins + 1, false);
509  }
510  else if (model_type=="ghosting") {
511  // set the view to either since it doesn't matter
512  // vertices==hyperedges
513  view_ = VERTEX_CENTRIC;
514  // unique set of global ids for the ghosts
515  typedef std::set<gno_t> ghost_t;
516 
517  // mapping from global id to the set of ghosts
518  typedef std::unordered_map<gno_t,ghost_t> ghost_map_t;
519 
520  primaryPinType=primaryEType;
521  adjacencyPinType =ia->getSecondAdjacencyEntityType();
522 
523  // number of layers of ghosting to do
524  unsigned int layers=2;
525  const Teuchos::ParameterEntry *pe3 = pl.getEntryPtr("ghost_layers");
526  if (pe3){
527  int l;
528  l = pe3->getValue<int>(&l);
529  layers = static_cast<unsigned int>(l);
530  }
531 
532  typedef int nonzero_t; // adjacency matrix doesn't need scalar_t
533  typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
534 
535  // Get an adjacency matrix representing the graph on the mesh
536  // using second adjacencies. If second adjacencies are not
537  // provided build the matrix from first adjacencies.
538  RCP<sparse_matrix_type> secondAdj;
539  if (!ia->avail2ndAdjs(primaryPinType,adjacencyPinType)) {
540  secondAdj=Zoltan2::get2ndAdjsMatFromAdjs<user_t>(ia,comm_,primaryPinType, adjacencyPinType);
541  }
542  else {
543  const offset_t* offsets;
544  const gno_t* adjacencyIds;
545  ia->get2ndAdjsView(primaryPinType,adjacencyPinType,offsets,adjacencyIds);
546  if (unique) {
547  Tpetra::global_size_t numGlobalCoords =
548  Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
549  oneToOneMap = rcp(new map_t(numGlobalCoords, gids_(), 0, comm));
550  // TODO KDD 1/17 It would be better to use minimum GID rather than
551  // TODO zero in the above Tpetra::Map constructor. Github #1024
552  }
553  secondAdj = rcp(new sparse_matrix_type(oneToOneMap,0));
554  for (size_t i=0; i<numLocalVertices_;i++) {
555  if (!isOwner_[i])
556  continue;
557  gno_t row = gids_[i];
558  offset_t num_adjs = offsets[i+1]-offsets[i];
559  ArrayRCP<nonzero_t> ones(num_adjs,1);
560  ArrayRCP<const gno_t> cols(adjacencyIds,offsets[i],num_adjs,false);
561  secondAdj->insertGlobalValues(row,cols(),ones());
562  }
563  secondAdj->fillComplete();
564  }
565 
566  //The mapping of the ghosts per hypergraph vertex
567  ghost_map_t ghosts;
568 
569  //Read the 1 layer ghosts from the second adjacency matrix
570  Array<gno_t> Indices;
571  Array<nonzero_t> Values;
572  for (unsigned int i=0;i<numLocalEdges_;i++) {
573  if (!isOwner_[i])
574  continue;
575  gno_t gid = edgeGids_[i];
576  size_t NumEntries = secondAdj->getNumEntriesInGlobalRow (gid);
577  Indices.resize (NumEntries);
578  Values.resize (NumEntries);
579  secondAdj->getGlobalRowCopy(gid,Indices(),Values(),NumEntries);
580  for (size_t j = 0; j < NumEntries; ++j) {
581  if(gid != Indices[j]) {
582  ghosts[gid].insert(Indices[j]);
583  }
584  }
585  }
586 
587  // The ith power of the second adjacency matrix is the ith layer of ghosts.
588  // Here we compute the ith power of the matrix and add the ith layer ghosts
589  // from the new matrix.
590  RCP<sparse_matrix_type> mat_old = secondAdj;
591  for (unsigned int i=1;i<layers;i++) {
592  RCP<sparse_matrix_type> mat_new =
593  rcp (new sparse_matrix_type(secondAdj->getRowMap(),0));
594  Tpetra::MatrixMatrix::Multiply(*mat_old,false,*secondAdj,false,*mat_new);
595  for (unsigned int j=0;j<numLocalEdges_;j++) {
596  if (!isOwner_[j])
597  continue;
598  gno_t gid = edgeGids_[j];
599  size_t NumEntries = mat_new->getNumEntriesInGlobalRow (gid);
600  Indices.resize(NumEntries);
601  Values.resize(NumEntries);
602  mat_new->getGlobalRowCopy(gid,Indices(),Values(),NumEntries);
603  for (size_t k = 0; k < NumEntries; ++k)
604  if(gid != Indices[k])
605  ghosts[gid].insert(Indices[k]);
606 
607  }
608  mat_old = mat_new;
609  }
610 
611  //Make the pins from the ghosts
612  for (size_t i=0;i<numLocalVertices_;i++) {//for each local entity
613  numLocalPins_+=ghosts[gids_[i]].size();
614  }
615  gno_t* temp_pins = new gno_t[numLocalPins_];
616  offset_t* temp_offsets = new offset_t[numLocalVertices_+1];
617  gno_t j=0;
618  for (size_t i=0;i<numLocalVertices_;i++) {//for each local entity
619  temp_offsets[i]=j;
620  if (!isOwner_[i])
621  continue;
622  typename ghost_t::iterator itr;
623  for (itr=ghosts[gids_[i]].begin();itr!=ghosts[gids_[i]].end();itr++) { //for each ghost of this entity
624  temp_pins[j]=*itr;
625  j++;
626 
627  }
628  }
629  temp_offsets[numLocalVertices_]=numLocalPins_;
630  pinGids_ = arcp<const gno_t>(temp_pins,0,numLocalPins_,true);
631  offsets_ = arcp<const offset_t>(temp_offsets,0,numLocalVertices_+1,true);
632 
633  //==============================Ghosting complete=================================
634  }
635 
636 
637  //Get the vertex weights
638  numWeightsPerVertex_ = ia->getNumWeightsPerID();
639 
640  if (numWeightsPerVertex_ > 0){
641  input_t *weightInfo = new input_t [numWeightsPerVertex_];
642  env_->localMemoryAssertion(__FILE__, __LINE__, numWeightsPerVertex_,
643  weightInfo);
644 
645  for (int idx=0; idx < numWeightsPerVertex_; idx++){
646  bool useNumNZ = ia->useDegreeAsWeight(idx);
647  if (useNumNZ){
648  scalar_t *wgts = new scalar_t [numLocalVertices_];
649  env_->localMemoryAssertion(__FILE__, __LINE__, numLocalVertices_, wgts);
650  ArrayRCP<const scalar_t> wgtArray =
651  arcp(wgts, 0, numLocalVertices_, true);
652  for (size_t i=0; i < numLocalVertices_; i++){
653  wgts[i] = offsets_[i+1] - offsets_[i];
654  }
655  weightInfo[idx] = input_t(wgtArray, 1);
656  }
657  else{
658  const scalar_t *weights=NULL;
659  int stride=0;
660  ia->getWeightsView(weights, stride, idx);
661  ArrayRCP<const scalar_t> wgtArray = arcp(weights, 0,
662  stride*numLocalVertices_,
663  false);
664  weightInfo[idx] = input_t(wgtArray, stride);
665  }
666  }
667 
668  vWeights_ = arcp<input_t>(weightInfo, 0, numWeightsPerVertex_, true);
669  }
670 
671  //TODO get the weights for edges, and pins(?)
672 
673  //Get the vertex coordinates from the primary types
674  typedef MeshAdapter<user_t> adapterWithCoords_t;
675  shared_GetVertexCoords<adapterWithCoords_t>(&(*ia));
676 
677  env_->timerStop(MACRO_TIMERS, "HyperGraphModel constructed from MeshAdapter");
678  print();
679 }
680 
682 
683 template <typename Adapter>
684 template <typename AdapterWithCoords>
685 void HyperGraphModel<Adapter>::shared_GetVertexCoords(const AdapterWithCoords *ia)
686 {
687  // get Vertex coordinates from input adapter
688 
689  vCoordDim_ = ia->getDimension();
690 
691  if (vCoordDim_ > 0){
692  input_t *coordInfo = new input_t [vCoordDim_];
693  env_->localMemoryAssertion(__FILE__, __LINE__, vCoordDim_, coordInfo);
694 
695  for (int dim=0; dim < vCoordDim_; dim++){
696  const scalar_t *coords=NULL;
697  int stride=0;
698  ia->getCoordinatesView(coords, stride, dim);
699  ArrayRCP<const scalar_t> coordArray = arcp(coords, 0,
700  stride*numLocalVertices_,
701  false);
702  coordInfo[dim] = input_t(coordArray, stride);
703  }
704 
705  vCoords_ = arcp<input_t>(coordInfo, 0, vCoordDim_, true);
706  }
707 }
708 
710 template <typename Adapter>
711 void HyperGraphModel<Adapter>::print()
712 {
713  //only prints the model if debug status is verbose
714  if (env_->getDebugLevel() < VERBOSE_DETAILED_STATUS)
715  return;
716 
717  std::ostream *os = env_->getDebugOStream();
718 
719  int me = comm_->getRank();
720  std::string fn(" ");
721 
722  *os << me << fn
723  << " Nvtx " << gids_.size()
724  << " Nedge " << edgeGids_.size()
725  << " NPins " << numLocalPins_
726  << " NVWgt " << numWeightsPerVertex_
727  << " NEWgt " << nWeightsPerEdge_
728  << " NPWgt " << nWeightsPerPin_
729  << " CDim " << vCoordDim_
730  << std::endl;
731 
732  for (lno_t i = 0; i < gids_.size(); i++) {
733  *os << me << fn << i << " VTXGID " << gids_[i]<<" isOwner: "<<isOwner_[i];
734  if (numWeightsPerVertex_==1)
735  *os << " weight: " << vWeights_[0][i];
736  if (view_==VERTEX_CENTRIC) {
737  *os <<" pins:";
738  for (offset_t j = offsets_[i]; j< offsets_[i+1];j++)
739  *os <<" "<<pinGids_[j];
740  }
741  *os << std::endl;
742  }
743  for (lno_t i = 0; i<edgeGids_.size();i++) {
744  *os << me << fn << i << " EDGEGID " << edgeGids_[i];
745  if (view_==HYPEREDGE_CENTRIC) {
746  *os <<":";
747  for (offset_t j = offsets_[i]; j< offsets_[i+1];j++)
748  *os <<" "<<pinGids_[j];
749  }
750  *os << std::endl;
751  }
752  if (vCoordDim_) {
753  for (lno_t i = 0; i < gids_.size(); i++) {
754  *os << me << fn << i << " COORDS " << gids_[i] << ": ";
755  for (int j = 0; j < vCoordDim_; j++)
756  *os << vCoords_[j][i] << " ";
757  *os << std::endl;
758  }
759  }
760  else
761  *os << me << fn << "NO COORDINATES AVAIL " << std::endl;
762 }
763 
764 } // namespace Zoltan2
765 
766 
767 #endif
768 
Zoltan2::MatrixAdapter
MatrixAdapter defines the adapter interface for matrices.
Definition: Zoltan2_MatrixAdapter.hpp:106
Z2_FORWARD_EXCEPTIONS
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Definition: Zoltan2_Exceptions.hpp:106
Zoltan2::VectorAdapter
VectorAdapter defines the interface for vector input.
Definition: Zoltan2_VectorAdapter.hpp:98
Zoltan2_InputTraits.hpp
Traits for application input objects.
Zoltan2::HyperGraphModel::HyperGraphModel
HyperGraphModel(const RCP< const VectorAdapter< userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags, CentricView view)
Definition: Zoltan2_HyperGraphModel.hpp:138
Zoltan2::HyperGraphModel::getVertexCoords
size_t getVertexCoords(ArrayView< input_t > &xyz) const
Sets pointers to this process' vertex coordinates, if available.
Definition: Zoltan2_HyperGraphModel.hpp:225
Zoltan2::HyperGraphModel::getLocalNumHyperEdges
size_t getLocalNumHyperEdges() const
Returns the number of hyper edges on this process. These are all hyper edges that have an adjacency t...
Definition: Zoltan2_HyperGraphModel.hpp:177
Zoltan2_ModelHelpers.hpp
Defines helper functions for use in the models.
Zoltan2::Model
The base class for all model classes.
Definition: Zoltan2_Model.hpp:110
Zoltan2::HyperGraphModel::getNumWeightsPerVertex
int getNumWeightsPerVertex() const
Returns the number (0 or greater) of weights per vertex.
Definition: Zoltan2_HyperGraphModel.hpp:189
Zoltan2::MeshEntityType
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
Definition: Zoltan2_MeshAdapter.hpp:63
Zoltan2::MACRO_TIMERS
Time an algorithm (or other entity) as a whole.
Definition: Zoltan2_Parameters.hpp:120
Zoltan2::HyperGraphModel::getGlobalNumObjects
size_t getGlobalNumObjects() const
Return the global number of objects.
Definition: Zoltan2_HyperGraphModel.hpp:305
Zoltan2_MatrixAdapter.hpp
Defines the MatrixAdapter interface.
Zoltan2::VERTEX_CENTRIC
Definition: Zoltan2_HyperGraphModel.hpp:76
Zoltan2::HyperGraphModel::getGlobalNumVertices
size_t getGlobalNumVertices() const
Returns the global number vertices.
Definition: Zoltan2_HyperGraphModel.hpp:171
Zoltan2::HyperGraphModel::getVertexMaps
void getVertexMaps(Teuchos::RCP< const map_t > &copiesMap, Teuchos::RCP< const map_t > &onetooneMap) const
Sets pointers to the vertex map with copies and the vertex map without copies Note: the pointers will...
Definition: Zoltan2_HyperGraphModel.hpp:252
Zoltan2::MeshAdapter
MeshAdapter defines the interface for mesh input.
Definition: Zoltan2_MeshAdapter.hpp:124
Zoltan2_VectorAdapter.hpp
Defines the VectorAdapter interface.
Zoltan2::HyperGraphModel::getLocalNumPins
size_t getLocalNumPins() const
Returns the local number of pins.
Definition: Zoltan2_HyperGraphModel.hpp:185
Zoltan2::HyperGraphModel::HyperGraphModel
HyperGraphModel(const RCP< const GraphAdapter< user_t, userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &modelFlags, CentricView view)
Definition: Zoltan2_HyperGraphModel.hpp:127
Zoltan2::HyperGraphModel::getEdgeList
size_t getEdgeList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
Sets pointers to this process' hyperedge Ids and their weights.
Definition: Zoltan2_HyperGraphModel.hpp:266
Zoltan2::VERBOSE_DETAILED_STATUS
include more detail about sub-steps
Definition: Zoltan2_Parameters.hpp:102
Zoltan2::HyperGraphModel::areVertexIDsUnique
bool areVertexIDsUnique() const
Returns true if the vertices are unique false otherwise.
Definition: Zoltan2_HyperGraphModel.hpp:159
Zoltan2::HyperGraphModel::getNumWeightsPerHyperEdge
int getNumWeightsPerHyperEdge() const
Returns the number (0 or greater) of weights per edge.
Definition: Zoltan2_HyperGraphModel.hpp:193
Zoltan2::HyperGraphModel::getNumWeightesPerPin
int getNumWeightesPerPin() const
Returns the number (0 or greater) of weights per pins.
Definition: Zoltan2_HyperGraphModel.hpp:197
Zoltan2::HyperGraphModel::getCentricView
CentricView getCentricView() const
Returns the centric view of the hypergraph.
Definition: Zoltan2_HyperGraphModel.hpp:155
Zoltan2::global_size_t
Tpetra::global_size_t global_size_t
Definition: Zoltan2_Standards.hpp:119
Zoltan2::HyperGraphModel::getLocalNumVertices
size_t getLocalNumVertices() const
Returns the number vertices on this process.
Definition: Zoltan2_HyperGraphModel.hpp:163
Zoltan2::HyperGraphModel::getLocalNumObjects
size_t getLocalNumObjects() const
Return the local number of objects.
Definition: Zoltan2_HyperGraphModel.hpp:303
Zoltan2_StridedData.hpp
This file defines the StridedData class.
Zoltan2::HyperGraphModel::getPinList
size_t getPinList(ArrayView< const gno_t > &pinIds, ArrayView< const offset_t > &offsets, ArrayView< input_t > &wgts) const
Sets pointers to this process' pins global Ids based on the centric view given by getCentricView()
Definition: Zoltan2_HyperGraphModel.hpp:288
Zoltan2::HyperGraphModel::getVertexList
size_t getVertexList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
Sets pointers to this process' vertex Ids and their weights.
Definition: Zoltan2_HyperGraphModel.hpp:210
Zoltan2::HyperGraphModel::getLocalNumOwnedVertices
size_t getLocalNumOwnedVertices() const
Returns the number vertices on this process that are owned.
Definition: Zoltan2_HyperGraphModel.hpp:167
Zoltan2_MeshAdapter.hpp
Defines the MeshAdapter interface.
Zoltan2::HyperGraphModel::HyperGraphModel
HyperGraphModel(const RCP< const MatrixAdapter< user_t, userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &modelFlags, CentricView view)
Constructor.
Definition: Zoltan2_HyperGraphModel.hpp:120
Zoltan2_GraphAdapter.hpp
Defines the GraphAdapter interface.
Zoltan2_Model.hpp
Defines the Model interface.
Zoltan2::HyperGraphModel::HyperGraphModel
HyperGraphModel(const RCP< const IdentifierAdapter< user_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags, CentricView view)
Definition: Zoltan2_HyperGraphModel.hpp:145
Zoltan2::HyperGraphModel::getOwnedList
size_t getOwnedList(ArrayView< bool > &isOwner) const
Sets pointer to the ownership of this processes vertices.
Definition: Zoltan2_HyperGraphModel.hpp:238
Zoltan2::HYPEREDGE_CENTRIC
Definition: Zoltan2_HyperGraphModel.hpp:75
weights
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Definition: rcbPerformanceZ1.cpp:82
Zoltan2
Definition: Zoltan2_AlgSerialGreedy.hpp:56
Zoltan2::HyperGraphModel
HyperGraphModel defines the interface required for hyper graph models.
Definition: Zoltan2_HyperGraphModel.hpp:90
Zoltan2::IdentifierAdapter
IdentifierAdapter defines the interface for identifiers.
Definition: Zoltan2_IdentifierAdapter.hpp:96
Zoltan2::GraphAdapter
GraphAdapter defines the interface for graph-based user data.
Definition: Zoltan2_GraphAdapter.hpp:99
Zoltan2_IdentifierAdapter.hpp
Defines the IdentifierAdapter interface.
Zoltan2::HyperGraphModel::getCoordinateDim
int getCoordinateDim() const
Returns the dimension (0 to 3) of vertex coordinates.
Definition: Zoltan2_HyperGraphModel.hpp:201
Zoltan2::modelFlag_t
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
Definition: Zoltan2_Model.hpp:93
Zoltan2::CentricView
CentricView
Enumerate the views for the pins: HYPEREDGE_CENTRIC: pins are the global ids of the vertices as seen ...
Definition: Zoltan2_HyperGraphModel.hpp:74
Zoltan2::StridedData
The StridedData class manages lists of weights or coordinates.
Definition: Zoltan2_StridedData.hpp:76
Zoltan2::HyperGraphModel::~HyperGraphModel
~HyperGraphModel()
Destructor.
Definition: Zoltan2_HyperGraphModel.hpp:107
Zoltan2::HyperGraphModel::getGlobalNumHyperEdges
size_t getGlobalNumHyperEdges() const
Returns the global number hyper edges.
Definition: Zoltan2_HyperGraphModel.hpp:181