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