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::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;
120 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
123 throw std::runtime_error(
"Building HyperGraphModel from MatrixAdapter not implemented yet");
127 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
130 throw std::runtime_error(
"Building HyperGraphModel from GraphAdapter not implemented yet");
134 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
138 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
141 throw std::runtime_error(
"cannot build HyperGraphModel from VectorAdapter");
145 const RCP<const Environment> &env,
const RCP<
const Comm<int> > &comm,
148 throw std::runtime_error(
"cannot build HyperGraphModel from IdentifierAdapter");
210 ArrayView<const gno_t> &Ids,
211 ArrayView<input_t> &wgts)
const 213 size_t nv = gids_.size();
215 wgts = vWeights_.view(0, numWeightsPerVertex_);
226 size_t nv = gids_.size();
227 xyz = vCoords_.view(0, vCoordDim_);
239 size_t nv = isOwner_.size();
240 isOwner = isOwner_(0, nv);
252 Teuchos::RCP<const map_t>& onetooneMap)
const 254 copiesMap = mapWithCopies;
255 onetooneMap = oneToOneMap;
266 ArrayView<const gno_t> &Ids,
267 ArrayView<input_t> &wgts)
const 269 size_t nv = edgeGids_.size();
270 Ids = edgeGids_(0, nv);
271 wgts = eWeights_.view(0, nWeightsPerEdge_);
288 ArrayView<const lno_t> &offsets,
289 ArrayView<input_t> &wgts)
const 291 pinIds = pinGids_(0, numLocalPins_);
292 offsets = offsets_.view(0, offsets_.size());
293 wgts = pWeights_.view(0, nWeightsPerPin_);
294 return pinGids_.size();
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;}
315 template <
typename AdapterWithCoords>
316 void shared_GetVertexCoords(
const AdapterWithCoords *ia);
319 const RCP<const Environment > env_;
320 const RCP<const Comm<int> > comm_;
324 ArrayRCP<const gno_t> gids_;
326 ArrayRCP<bool> isOwner_;
328 int numWeightsPerVertex_;
329 ArrayRCP<input_t> vWeights_;
332 ArrayRCP<input_t> vCoords_;
334 ArrayRCP<const gno_t> edgeGids_;
336 int nWeightsPerEdge_;
337 ArrayRCP<input_t> eWeights_;
339 ArrayRCP<const gno_t> pinGids_;
340 ArrayRCP<const lno_t> offsets_;
343 ArrayRCP<input_t> pWeights_;
347 size_t numLocalVertices_;
348 size_t numOwnedVertices_;
349 size_t numGlobalVertices_;
350 size_t numLocalEdges_;
351 size_t numGlobalEdges_;
352 size_t numLocalPins_;
355 Teuchos::RCP<const map_t> mapWithCopies;
356 Teuchos::RCP<const map_t> oneToOneMap;
369 template <
typename Adapter>
372 const RCP<const Environment> &env,
373 const RCP<
const Comm<int> > &comm,
381 numWeightsPerVertex_(0),
392 numLocalVertices_(0),
393 numGlobalVertices_(0),
398 env_->timerStart(
MACRO_TIMERS,
"HyperGraphModel constructed from MeshAdapter");
408 std::string model_type(
"traditional");
409 const Teuchos::ParameterList &pl = env->getParameters();
410 const Teuchos::ParameterEntry *pe2 = pl.getEntryPtr(
"hypergraph_model_type");
412 model_type = pe2->getValue<std::string>(&model_type);
420 gno_t
const *vtxIds=NULL;
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_);
432 gids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_,
false);
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;
444 unique = ia->areEntityIDsUnique(ia->getPrimaryEntityType());
445 numOwnedVertices_=numLocalVertices_;
446 isOwner_ = ArrayRCP<bool>(numLocalVertices_,
true);
450 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
451 mapWithCopies = rcp(
new map_t(numGlobalCoords, gids_(), 0, comm));
454 oneToOneMap = Tpetra::createOneToOne<lno_t, gno_t>(mapWithCopies);
456 numOwnedVertices_=oneToOneMap->getNodeNumElements();
457 for (
size_t i=0;i<numLocalVertices_;i++) {
458 isOwner_[i] = oneToOneMap->isNodeGlobalElement(gids_[i]);
463 if (model_type==
"traditional") {
467 gno_t
const *edgeIds=NULL;
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_);
476 edgeGids_ = arcp<const gno_t>(edgeIds, 0, numLocalEdges_,
false);
478 else if (model_type==
"ghosting") {
480 numLocalEdges_ = numLocalVertices_;
481 edgeGids_ = arcp<const gno_t>(vtxIds, 0, numLocalVertices_,
false);
482 numGlobalEdges_ = numGlobalVertices_;
488 size_t numPrimaryPins = numLocalVertices_;
490 primaryPinType = adjacencyEType;
491 adjacencyPinType = primaryEType;
492 numPrimaryPins = numLocalEdges_;
494 if (model_type==
"traditional") {
496 gno_t
const *nborIds=NULL;
497 lno_t
const *offsets=NULL;
500 ia->getAdjsView(primaryPinType,adjacencyPinType,offsets,nborIds);
504 numLocalPins_ = offsets[numPrimaryPins];
506 pinGids_ = arcp<const gno_t>(nborIds, 0, numLocalPins_,
false);
507 offsets_ = arcp<const lno_t>(offsets, 0, numPrimaryPins + 1,
false);
509 else if (model_type==
"ghosting") {
514 typedef std::set<gno_t> ghost_t;
517 typedef std::unordered_map<gno_t,ghost_t> ghost_map_t;
519 primaryPinType=primaryEType;
520 adjacencyPinType =ia->getSecondAdjacencyEntityType();
523 unsigned int layers=2;
524 const Teuchos::ParameterEntry *pe3 = pl.getEntryPtr(
"ghost_layers");
527 l = pe3->getValue<
int>(&l);
528 layers =
static_cast<unsigned int>(l);
531 typedef int nonzero_t;
532 typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
537 RCP<sparse_matrix_type> secondAdj;
538 if (!ia->avail2ndAdjs(primaryPinType,adjacencyPinType)) {
539 secondAdj=Zoltan2::get2ndAdjsMatFromAdjs<user_t>(ia,comm_,primaryPinType, adjacencyPinType);
542 const lno_t* offsets;
543 const gno_t* adjacencyIds;
544 ia->get2ndAdjsView(primaryPinType,adjacencyPinType,offsets,adjacencyIds);
547 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
548 oneToOneMap = rcp(
new map_t(numGlobalCoords, gids_(), 0, comm));
552 secondAdj = rcp(
new sparse_matrix_type(oneToOneMap,0));
553 for (
size_t i=0; i<numLocalVertices_;i++) {
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());
562 secondAdj->fillComplete();
569 Array<gno_t> Indices;
570 Array<nonzero_t> Values;
571 for (
unsigned int i=0;i<numLocalEdges_;i++) {
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]);
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++) {
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]);
611 for (
size_t i=0;i<numLocalVertices_;i++) {
612 numLocalPins_+=ghosts[gids_[i]].size();
614 gno_t* temp_pins =
new gno_t[numLocalPins_];
615 lno_t* temp_offsets =
new lno_t[numLocalVertices_+1];
617 for (
size_t i=0;i<numLocalVertices_;i++) {
621 typename ghost_t::iterator itr;
622 for (itr=ghosts[gids_[i]].begin();itr!=ghosts[gids_[i]].end();itr++) {
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);
637 numWeightsPerVertex_ = ia->getNumWeightsPerID();
639 if (numWeightsPerVertex_ > 0){
640 input_t *weightInfo =
new input_t [numWeightsPerVertex_];
641 env_->localMemoryAssertion(__FILE__, __LINE__, numWeightsPerVertex_,
644 for (
int idx=0; idx < numWeightsPerVertex_; idx++){
645 bool useNumNZ = ia->useDegreeAsWeight(idx);
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];
654 weightInfo[idx] = input_t(wgtArray, 1);
659 ia->getWeightsView(weights, stride, idx);
660 ArrayRCP<const scalar_t> wgtArray = arcp(weights, 0,
661 stride*numLocalVertices_,
663 weightInfo[idx] = input_t(wgtArray, stride);
667 vWeights_ = arcp<input_t>(weightInfo, 0, numWeightsPerVertex_,
true);
674 shared_GetVertexCoords<adapterWithCoords_t>(&(*ia));
676 env_->timerStop(
MACRO_TIMERS,
"HyperGraphModel constructed from MeshAdapter");
682 template <
typename Adapter>
683 template <
typename AdapterWithCoords>
688 vCoordDim_ = ia->getDimension();
691 input_t *coordInfo =
new input_t [vCoordDim_];
692 env_->localMemoryAssertion(__FILE__, __LINE__, vCoordDim_, coordInfo);
694 for (
int dim=0; dim < vCoordDim_; dim++){
695 const scalar_t *coords=NULL;
697 ia->getCoordinatesView(coords, stride, dim);
698 ArrayRCP<const scalar_t> coordArray = arcp(coords, 0,
699 stride*numLocalVertices_,
701 coordInfo[dim] = input_t(coordArray, stride);
704 vCoords_ = arcp<input_t>(coordInfo, 0, vCoordDim_,
true);
709 template <
typename Adapter>
716 std::ostream *os = env_->getDebugOStream();
718 int me = comm_->getRank();
722 <<
" Nvtx " << gids_.size()
723 <<
" Nedge " << edgeGids_.size()
724 <<
" NPins " << numLocalPins_
725 <<
" NVWgt " << numWeightsPerVertex_
726 <<
" NEWgt " << nWeightsPerEdge_
727 <<
" NPWgt " << nWeightsPerPin_
728 <<
" CDim " << vCoordDim_
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];
737 for (lno_t j = offsets_[i]; j< offsets_[i+1];j++)
738 *os <<
" "<<pinGids_[j];
742 for (lno_t i = 0; i<edgeGids_.size();i++) {
743 *os << me << fn << i <<
" EDGEGID " << edgeGids_[i];
746 for (lno_t j = offsets_[i]; j< offsets_[i+1];j++)
747 *os <<
" "<<pinGids_[j];
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] <<
" ";
760 *os << me << fn <<
"NO COORDINATES AVAIL " << std::endl;
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' 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.
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' 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.
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' 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...
~HyperGraphModel()
Destructor.
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' vertex coordinates, if available.