50 #ifndef _ZOLTAN2_HYPERGRAPHMODEL_HPP_
51 #define _ZOLTAN2_HYPERGRAPHMODEL_HPP_
64 #include <unordered_map>
66 #include <Teuchos_Hashtable.hpp>
89 template <
typename Adapter>
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;
121 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
124 throw std::runtime_error(
"Building HyperGraphModel from MatrixAdapter not implemented yet");
128 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
131 throw std::runtime_error(
"Building HyperGraphModel from GraphAdapter not implemented yet");
135 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
139 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
142 throw std::runtime_error(
"cannot build HyperGraphModel from VectorAdapter");
146 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
149 throw std::runtime_error(
"cannot build HyperGraphModel from IdentifierAdapter");
211 ArrayView<const gno_t> &Ids,
212 ArrayView<input_t> &wgts)
const
214 size_t nv = gids_.size();
216 wgts = vWeights_.view(0, numWeightsPerVertex_);
227 size_t nv = gids_.size();
228 xyz = vCoords_.view(0, vCoordDim_);
240 size_t nv = isOwner_.size();
241 isOwner = isOwner_(0, nv);
253 Teuchos::RCP<const map_t>& onetooneMap)
const
255 copiesMap = mapWithCopies;
256 onetooneMap = oneToOneMap;
267 ArrayView<const gno_t> &Ids,
268 ArrayView<input_t> &wgts)
const
270 size_t nv = edgeGids_.size();
271 Ids = edgeGids_(0, nv);
272 wgts = eWeights_.view(0, nWeightsPerEdge_);
289 ArrayView<const offset_t> &offsets,
290 ArrayView<input_t> &wgts)
const
292 pinIds = pinGids_(0, numLocalPins_);
293 offsets = offsets_.view(0, offsets_.size());
294 wgts = pWeights_.view(0, nWeightsPerPin_);
295 return pinGids_.size();
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;}
316 template <
typename AdapterWithCoords>
317 void shared_GetVertexCoords(
const AdapterWithCoords *ia);
320 const RCP<const Environment > env_;
321 const RCP<const Comm<int> > comm_;
325 ArrayRCP<const gno_t> gids_;
327 ArrayRCP<bool> isOwner_;
329 int numWeightsPerVertex_;
330 ArrayRCP<input_t> vWeights_;
333 ArrayRCP<input_t> vCoords_;
335 ArrayRCP<const gno_t> edgeGids_;
337 int nWeightsPerEdge_;
338 ArrayRCP<input_t> eWeights_;
340 ArrayRCP<const gno_t> pinGids_;
341 ArrayRCP<const offset_t> offsets_;
344 ArrayRCP<input_t> pWeights_;
348 size_t numLocalVertices_;
349 size_t numOwnedVertices_;
350 size_t numGlobalVertices_;
351 size_t numLocalEdges_;
352 size_t numGlobalEdges_;
353 size_t numLocalPins_;
356 Teuchos::RCP<const map_t> mapWithCopies;
357 Teuchos::RCP<const map_t> oneToOneMap;
370 template <
typename Adapter>
373 const RCP<const Environment> &env,
374 const RCP<
const Comm<int> > &comm,
382 numWeightsPerVertex_(0),
393 numLocalVertices_(0),
394 numGlobalVertices_(0),
399 env_->timerStart(
MACRO_TIMERS,
"HyperGraphModel constructed from MeshAdapter");
409 std::string model_type(
"traditional");
410 const Teuchos::ParameterList &pl = env->getParameters();
411 const Teuchos::ParameterEntry *pe2 = pl.getEntryPtr(
"hypergraph_model_type");
413 model_type = pe2->getValue<std::string>(&model_type);
421 gno_t
const *vtxIds=NULL;
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_);
433 gids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_,
false);
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;
445 unique = ia->areEntityIDsUnique(ia->getPrimaryEntityType());
446 numOwnedVertices_=numLocalVertices_;
447 isOwner_ = ArrayRCP<bool>(numLocalVertices_,
true);
451 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
452 mapWithCopies = rcp(
new map_t(numGlobalCoords, gids_(), 0, comm));
455 oneToOneMap = Tpetra::createOneToOne<lno_t, gno_t>(mapWithCopies);
457 numOwnedVertices_=oneToOneMap->getNodeNumElements();
458 for (
size_t i=0;i<numLocalVertices_;i++) {
459 isOwner_[i] = oneToOneMap->isNodeGlobalElement(gids_[i]);
464 if (model_type==
"traditional") {
468 gno_t
const *edgeIds=NULL;
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_);
477 edgeGids_ = arcp<const gno_t>(edgeIds, 0, numLocalEdges_,
false);
479 else if (model_type==
"ghosting") {
481 numLocalEdges_ = numLocalVertices_;
482 edgeGids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_,
false);
483 numGlobalEdges_ = numGlobalVertices_;
489 size_t numPrimaryPins = numLocalVertices_;
491 primaryPinType = adjacencyEType;
492 adjacencyPinType = primaryEType;
493 numPrimaryPins = numLocalEdges_;
495 if (model_type==
"traditional") {
497 gno_t
const *nborIds=NULL;
498 offset_t
const *offsets=NULL;
501 ia->getAdjsView(primaryPinType,adjacencyPinType,offsets,nborIds);
505 numLocalPins_ = offsets[numPrimaryPins];
507 pinGids_ = arcp<const gno_t>(nborIds, 0, numLocalPins_,
false);
508 offsets_ = arcp<const offset_t>(offsets, 0, numPrimaryPins + 1,
false);
510 else if (model_type==
"ghosting") {
515 typedef std::set<gno_t> ghost_t;
518 typedef std::unordered_map<gno_t,ghost_t> ghost_map_t;
520 primaryPinType=primaryEType;
521 adjacencyPinType =ia->getSecondAdjacencyEntityType();
524 unsigned int layers=2;
525 const Teuchos::ParameterEntry *pe3 = pl.getEntryPtr(
"ghost_layers");
528 l = pe3->getValue<
int>(&l);
529 layers = static_cast<unsigned int>(l);
532 typedef int nonzero_t;
533 typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
538 RCP<sparse_matrix_type> secondAdj;
539 if (!ia->avail2ndAdjs(primaryPinType,adjacencyPinType)) {
540 secondAdj=Zoltan2::get2ndAdjsMatFromAdjs<user_t>(ia,comm_,primaryPinType, adjacencyPinType);
543 const offset_t* offsets;
544 const gno_t* adjacencyIds;
545 ia->get2ndAdjsView(primaryPinType,adjacencyPinType,offsets,adjacencyIds);
548 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
549 oneToOneMap = rcp(
new map_t(numGlobalCoords, gids_(), 0, comm));
553 secondAdj = rcp(
new sparse_matrix_type(oneToOneMap,0));
554 for (
size_t i=0; i<numLocalVertices_;i++) {
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());
563 secondAdj->fillComplete();
570 Array<gno_t> Indices;
571 Array<nonzero_t> Values;
572 for (
unsigned int i=0;i<numLocalEdges_;i++) {
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]);
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++) {
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]);
612 for (
size_t i=0;i<numLocalVertices_;i++) {
613 numLocalPins_+=ghosts[gids_[i]].size();
615 gno_t* temp_pins =
new gno_t[numLocalPins_];
616 offset_t* temp_offsets =
new offset_t[numLocalVertices_+1];
618 for (
size_t i=0;i<numLocalVertices_;i++) {
622 typename ghost_t::iterator itr;
623 for (itr=ghosts[gids_[i]].begin();itr!=ghosts[gids_[i]].end();itr++) {
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);
638 numWeightsPerVertex_ = ia->getNumWeightsPerID();
640 if (numWeightsPerVertex_ > 0){
641 input_t *weightInfo =
new input_t [numWeightsPerVertex_];
642 env_->localMemoryAssertion(__FILE__, __LINE__, numWeightsPerVertex_,
645 for (
int idx=0; idx < numWeightsPerVertex_; idx++){
646 bool useNumNZ = ia->useDegreeAsWeight(idx);
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];
655 weightInfo[idx] = input_t(wgtArray, 1);
660 ia->getWeightsView(
weights, stride, idx);
661 ArrayRCP<const scalar_t> wgtArray = arcp(
weights, 0,
662 stride*numLocalVertices_,
664 weightInfo[idx] = input_t(wgtArray, stride);
668 vWeights_ = arcp<input_t>(weightInfo, 0, numWeightsPerVertex_,
true);
675 shared_GetVertexCoords<adapterWithCoords_t>(&(*ia));
677 env_->timerStop(
MACRO_TIMERS,
"HyperGraphModel constructed from MeshAdapter");
683 template <
typename Adapter>
684 template <
typename AdapterWithCoords>
689 vCoordDim_ = ia->getDimension();
692 input_t *coordInfo =
new input_t [vCoordDim_];
693 env_->localMemoryAssertion(__FILE__, __LINE__, vCoordDim_, coordInfo);
695 for (
int dim=0; dim < vCoordDim_; dim++){
696 const scalar_t *coords=NULL;
698 ia->getCoordinatesView(coords, stride, dim);
699 ArrayRCP<const scalar_t> coordArray = arcp(coords, 0,
700 stride*numLocalVertices_,
702 coordInfo[dim] = input_t(coordArray, stride);
705 vCoords_ = arcp<input_t>(coordInfo, 0, vCoordDim_,
true);
710 template <
typename Adapter>
711 void HyperGraphModel<Adapter>::print()
717 std::ostream *os = env_->getDebugOStream();
719 int me = comm_->getRank();
723 <<
" Nvtx " << gids_.size()
724 <<
" Nedge " << edgeGids_.size()
725 <<
" NPins " << numLocalPins_
726 <<
" NVWgt " << numWeightsPerVertex_
727 <<
" NEWgt " << nWeightsPerEdge_
728 <<
" NPWgt " << nWeightsPerPin_
729 <<
" CDim " << vCoordDim_
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];
738 for (offset_t j = offsets_[i]; j< offsets_[i+1];j++)
739 *os <<
" "<<pinGids_[j];
743 for (lno_t i = 0; i<edgeGids_.size();i++) {
744 *os << me << fn << i <<
" EDGEGID " << edgeGids_[i];
747 for (offset_t j = offsets_[i]; j< offsets_[i+1];j++)
748 *os <<
" "<<pinGids_[j];
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] <<
" ";
761 *os << me << fn <<
"NO COORDINATES AVAIL " << std::endl;