Isorropia: Partitioning, Load Balancing and more
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Isorropia::Epetra::CostDescriber Class Reference

#include <Isorropia_EpetraCostDescriber.hpp>

Inheritance diagram for Isorropia::Epetra::CostDescriber:
Inheritance graph
[legend]
Collaboration diagram for Isorropia::Epetra::CostDescriber:
Collaboration graph
[legend]

Public Member Functions

 CostDescriber ()
 Constructor. More...
 
 ~CostDescriber ()
 Destructor. More...
 
 CostDescriber (const CostDescriber &costs)
 Copy constructor. More...
 
void setVertexWeights (Teuchos::RCP< const Epetra_Vector > vwgts)
 setVertexWeights is called by a process to supply the weight of each vertex (row) in the graph or hypergraph. More...
 
void setGraphEdgeWeights (Teuchos::RCP< const Epetra_CrsMatrix > gewts)
 setGraphEdgeWeights is called by a process to supply the weights for each of the edges of its vertices. More...
 
void setHypergraphEdgeWeights (Teuchos::RCP< const Epetra_Vector > hgewts)
 setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix. More...
 
void setHypergraphEdgeWeights (int numHGedges, const int *hgGIDs, const float *hgEwgts)
 setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix. More...
 
void setHypergraphEdgeWeights (int numHGedges, const int *hgGIDs, const double *hgEwgts)
 setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix. More...
 
void getCosts (std::map< int, float > &vertexWeights, std::map< int, std::map< int, float > > &graphEdgeWeights, std::map< int, float > &hypergraphEdgeWeights) const
 Get the contents of this CostDescriber. More...
 
void show_cd (std::ostream &) const
 Print out the contents of this CostDescriber. More...
 

Private Member Functions

void _transformWeights (const Epetra_Import &importer)
 
int _compareBeforeAndAfterGraph (const Epetra_RowMatrix *in_m, const Epetra_RowMatrix *out_m, const Epetra_CrsGraph *in_g, const Epetra_CrsGraph *out_g, const Epetra_Import &importer, std::vector< double > &balance, std::vector< int > &numCuts, std::vector< double > &cutWgt, std::vector< double > &cutn, std::vector< double > &cutl) const
 
void setParameters (const Teuchos::ParameterList &paramlist)
 Set parameters for the CostDescriber instance. More...
 
bool haveVertexWeights () const
 Query whether vertex weights have been supplied by the application. More...
 
int getNumVertices () const
 Get the number of vertices for which this process supplied vertex weights. More...
 
void getVertexWeights (int numVertices, int *global_ids, float *weights) const
 Get lists of the vertex ids and weights supplied by this process. More...
 
bool haveGraphEdgeWeights () const
 Query whether graph edge weights have been supplied by the application. More...
 
int getNumGraphEdges (int vertex_global_id) const
 Get the number of graph edges for a specified vertex. More...
 
int getGraphEdgeVertices (std::set< int > &gids) const
 Get the set of global IDs for the vertices that we have edge information for. More...
 
void getGraphEdgeWeights (int vertex_global_id, int num_neighbors, int *neighbor_global_ids, float *weights) const
 Get the graph edge weights for a specified vertex. More...
 
bool haveHypergraphEdgeWeights () const
 Query whether hypergraph edge weights have been supplied by the application. More...
 
int getNumHypergraphEdgeWeights () const
 Get the number of Hypergraph edges. More...
 
void getHypergraphEdgeWeights (int numEdges, int *global_ids, float *weights) const
 Get the hypergraph edge weights that were supplied by this process. More...
 
int getHypergraphEdgeWeights (std::map< int, float > &wgtMap) const
 Return the CostDescribers hypergraph edge weights as a map from hyperedge (column) global ID to weight. More...
 
int getVertexWeights (std::map< int, float > &wgtMap) const
 Get vertex weights in the form of a map from vertex global ID to vertex weight. More...
 
const Epetra_Vector & getVertexWeights ()
 
int getGraphEdgeWeights (int vertex_global_id, std::map< int, float > &wgtMap) const
 getGraphEdgeWeights is called to obtain the graph edge weights for a given vertex (row). More...
 
bool haveGlobalVertexWeights () const
 haveGlobalVertexWeights returns true if any process in the application has supplied vertex weights, it returns false otherwise. More...
 
void setNumGlobalVertexWeights (int num)
 setNumGlobalVertexWeights may be used to set the count of the global number of vertex weights supplied to the CostDescriber More...
 
bool haveGlobalGraphEdgeWeights () const
 haveGlobalGraphEdgeWeights returns true if any process in the application has supplied graph edge weights, it returns false otherwise. More...
 
void setNumGlobalGraphEdgeWeights (int num)
 setNumGlobalGraphEdgeWeights may be used to set the count of the global number of graph edge weights supplied to the CostDescriber More...
 
bool haveGlobalHypergraphEdgeWeights () const
 haveGlobalHypergraphEdgeWeights returns true if any process in the application has supplied hyperedge weights, it returns false otherwise. More...
 
void setNumGlobalHypergraphEdgeWeights (int num)
 setNumGlobalHypergraphEdgeWeights may be used to set the count of the global number of hyperedge weights supplied to the CostDescriber More...
 
void allocate_hg_edge_weights_ (int n)
 Dynamically allocate storage for hypergraph edge weights. More...
 
void free_hg_edge_weights_ ()
 Free storage used by hypergraph edge weights. More...
 
int getEdges (int vertexGID, int len, int *nborGID, float *weights) const
 getEdges creates an array of the neighbors and edge weights for given vertex. More...
 

Private Attributes

Teuchos::RCP< const Epetra_Vector > vertex_weights_
 
Teuchos::RCP< const
Epetra_CrsMatrix > 
graph_edge_weights_
 
std::set< int > graph_self_edges_
 
Teuchos::ParameterList paramlist_
 
int * hg_edge_gids_
 
float * hg_edge_weights_
 
int num_hg_edge_weights_
 
int numGlobalVertexWeights_
 
int numGlobalGraphEdgeWeights_
 
int numGlobalHypergraphEdgeWeights_
 

Friends

class Isorropia::Operator
 
class Isorropia::Epetra::ZoltanLib::QueryObject
 
class Isorropia::Epetra::ZoltanLibClass
 
std::ostream & operator<< (std::ostream &, const Isorropia::Epetra::CostDescriber &cd)
 Overloaded << operator for CostDescriber object. More...
 

Detailed Description

Examples:
geometric/example_rcb.cpp, graphedge_weights.cpp, hgedge_weights.cpp, matrix_1.cpp, part_redist.cpp, and vert_weights.cpp.

Constructor & Destructor Documentation

Isorropia::Epetra::CostDescriber::CostDescriber ( )

Constructor.

Isorropia::Epetra::CostDescriber::~CostDescriber ( )
virtual

Destructor.

Reimplemented from Isorropia::CostDescriber.

Isorropia::Epetra::CostDescriber::CostDescriber ( const CostDescriber costs)

Copy constructor.

Member Function Documentation

void Isorropia::Epetra::CostDescriber::setVertexWeights ( Teuchos::RCP< const Epetra_Vector >  vwgts)

setVertexWeights is called by a process to supply the weight of each vertex (row) in the graph or hypergraph.

If the object to be partitioned is instead an Epetra_MultiVector representing real coordinates, then the weights represent the weight assigned to each coordinate.

Parameters
vwgts[in] vector of weights, one for each vertex
Examples:
geometric/example_rcb.cpp, and vert_weights.cpp.
void Isorropia::Epetra::CostDescriber::setGraphEdgeWeights ( Teuchos::RCP< const Epetra_CrsMatrix >  gewts)

setGraphEdgeWeights is called by a process to supply the weights for each of the edges of its vertices.

An edge corresponds to a non-zero in the row representing the vertex.

This method is called only when performing graph partitioning with a square symmetric matrix. For hypergraph partitioning call the equivalent hypergraph method.

Parameters
gewtsan Epetra_CrsMatrix supplied by the application, each non-zero represents a weight for an edge
Examples:
graphedge_weights.cpp.
void Isorropia::Epetra::CostDescriber::setHypergraphEdgeWeights ( Teuchos::RCP< const Epetra_Vector >  hgewts)

setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix.

(A hyperedge can in general link more than one vertex.)

Matrices that represent hypergraphs are not in general square. There may be more or fewer hyperedges (columns) than vertices (rows).

Parameters
hgewtsan Epetra_Vector containing the weights for each hyperedge.
Examples:
hgedge_weights.cpp.
void Isorropia::Epetra::CostDescriber::setHypergraphEdgeWeights ( int  numHGedges,
const int *  hgGIDs,
const float *  hgEwgts 
)

setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix.

(A hyperedge can in general link more than one vertex.)

Matrices that represent hypergraphs are not in general square. There may be more or fewer hyperedges (columns) than vertices (rows).

More than one process can supply a weight for the same column. (So there is no concept of a process owning a hyperedge.) Zoltan combines these weights according to the setting of the PHG_EDGE_WEIGHT_OPERATION parameter.

Parameters
numHGedgesthe length of the hgGIDs and heEwgts arrays
hgGIDsthe global ID for each hyperedge this process will supply a weight for
hgEwgtsthe hyperedge weight corresponding to each hyperedge listed in hgGIDs
void Isorropia::Epetra::CostDescriber::setHypergraphEdgeWeights ( int  numHGedges,
const int *  hgGIDs,
const double *  hgEwgts 
)

setHypergraphEdgeWeights is called by processes in an application to supply weights for the hyperedges, which are represented by the columns of the matrix.

(A hyperedge can in general link more than one vertex.)

Matrices that represent hypergraphs are not in general square. There may be more or fewer hyperedges (columns) than vertices (rows).

More than one process can supply a weight for the same column. (So there is no concept of a process owning a hyperedge.) Zoltan combines these weights according to the setting of the PHG_EDGE_WEIGHT_OPERATION parameter.

Parameters
numHGedgesthe length of the hgGIDs and heEwgts arrays
hgGIDsthe global ID for each hyperedge this process will supply a weight for
hgEwgtsthe hyperedge weight corresponding to each hyperedge listed in hgGIDs
void Isorropia::Epetra::CostDescriber::getCosts ( std::map< int, float > &  vertexWeights,
std::map< int, std::map< int, float > > &  graphEdgeWeights,
std::map< int, float > &  hypergraphEdgeWeights 
) const

Get the contents of this CostDescriber.

Parameters
vertexWeightsis set to a mapping from vertex global IDs to their weights
graphEdgeWeightsis a mapping from vertex global IDs to a map from neighboring IDs to edge weights
hypergraphEdgeWeightsis a mapping from hyperedge (column) global IDs to hyperedge weights
void Isorropia::Epetra::CostDescriber::show_cd ( std::ostream &  ) const

Print out the contents of this CostDescriber.

void Isorropia::Epetra::CostDescriber::_transformWeights ( const Epetra_Import &  importer)
private
int Isorropia::Epetra::CostDescriber::_compareBeforeAndAfterGraph ( const Epetra_RowMatrix *  in_m,
const Epetra_RowMatrix *  out_m,
const Epetra_CrsGraph *  in_g,
const Epetra_CrsGraph *  out_g,
const Epetra_Import &  importer,
std::vector< double > &  balance,
std::vector< int > &  numCuts,
std::vector< double > &  cutWgt,
std::vector< double > &  cutn,
std::vector< double > &  cutl 
) const
private
void Isorropia::Epetra::CostDescriber::setParameters ( const Teuchos::ParameterList &  paramlist)
privatevirtual

Set parameters for the CostDescriber instance.

The contents of the input paramlist object are copied into an internal ParameterList attribute. Instances of this interface should not retain a reference to the input ParameterList after this method returns.

Implements Isorropia::CostDescriber.

bool Isorropia::Epetra::CostDescriber::haveVertexWeights ( ) const
privatevirtual

Query whether vertex weights have been supplied by the application.

Returns
returns true if the application has supplied vertex weights with the CostDescriber, false otherwise

Implements Isorropia::CostDescriber.

int Isorropia::Epetra::CostDescriber::getNumVertices ( ) const
privatevirtual

Get the number of vertices for which this process supplied vertex weights.

Vertices typically correspond to matrix rows.

Returns
returns the number of vertices on this process for which the CostDescriber has vertex weights

Implements Isorropia::CostDescriber.

void Isorropia::Epetra::CostDescriber::getVertexWeights ( int  numVertices,
int *  global_ids,
float *  weights 
) const
privatevirtual

Get lists of the vertex ids and weights supplied by this process.

Parameters
numVerticessize of global_ids and weights arrays
global_idspointer to an array of vertex global IDs, allocated by the caller.
weightspointer to an array of vertex weights, allocated by the caller.

Implements Isorropia::CostDescriber.

bool Isorropia::Epetra::CostDescriber::haveGraphEdgeWeights ( ) const
privatevirtual

Query whether graph edge weights have been supplied by the application.

Returns
returns true if the application has supplied graph edge weights with the CostDescriber, false otherwise

Implements Isorropia::CostDescriber.

int Isorropia::Epetra::CostDescriber::getNumGraphEdges ( int  vertex_global_id) const
privatevirtual

Get the number of graph edges for a specified vertex.

Graph edges typically correspond to matrix nonzeros.

Parameters
vertex_global_idThe global ID for the vertex (on this process) for which the number of edges is desired
Returns
the number of graph edges supplied by this process for this vertex

Implements Isorropia::CostDescriber.

int Isorropia::Epetra::CostDescriber::getGraphEdgeVertices ( std::set< int > &  gids) const
private

Get the set of global IDs for the vertices that we have edge information for.

Parameters
gidswill be set to the global IDs of the vertices for which neighbor and edge weight information have been provided
void Isorropia::Epetra::CostDescriber::getGraphEdgeWeights ( int  vertex_global_id,
int  num_neighbors,
int *  neighbor_global_ids,
float *  weights 
) const
privatevirtual

Get the graph edge weights for a specified vertex.

Parameters
vertex_global_idvertex global ID (on this process) for which edge information is requested
num_neighborssize for which neighbor_global_ids and weights had been preallocated
neighbor_global_idsbuffer allocated by caller, on return will contain a list of neighbor vertex global IDs
weightsbuffer allocated by caller, on return will contain a weight for each edge indicated in neighbor_global_ids

Implements Isorropia::CostDescriber.

bool Isorropia::Epetra::CostDescriber::haveHypergraphEdgeWeights ( ) const
privatevirtual

Query whether hypergraph edge weights have been supplied by the application.

Returns
returns true if the application has supplied hypergraph edge weights with the CostDescriber, false otherwise

Implements Isorropia::CostDescriber.

int Isorropia::Epetra::CostDescriber::getNumHypergraphEdgeWeights ( ) const
privatevirtual

Get the number of Hypergraph edges.

Hypergraph edges typically correspond to matrix columns.

Returns
returns the number of hypergraph edge weights supplied by this process

Implements Isorropia::CostDescriber.

void Isorropia::Epetra::CostDescriber::getHypergraphEdgeWeights ( int  numEdges,
int *  global_ids,
float *  weights 
) const
privatevirtual

Get the hypergraph edge weights that were supplied by this process.

Parameters
numEdgessize for which global_ids and weights had been preallocated
global_idsbuffer allocated by caller, on return will contain a list of hyperedge global IDs
weightsbuffer allocated by caller, on return will contain a weight for each hyperedge indicated in global_ids

Implements Isorropia::CostDescriber.

int Isorropia::Epetra::CostDescriber::getHypergraphEdgeWeights ( std::map< int, float > &  wgtMap) const
private

Return the CostDescribers hypergraph edge weights as a map from hyperedge (column) global ID to weight.

Parameters
wgtMapwill be set to a map from hyperedge global ID to hyperedge weight
int Isorropia::Epetra::CostDescriber::getVertexWeights ( std::map< int, float > &  wgtMap) const
private

Get vertex weights in the form of a map from vertex global ID to vertex weight.

Parameters
wgtMapa map supplied by the caller, the vertex weights will be added to this map
Returns
the size of wgtMap
const Epetra_Vector& Isorropia::Epetra::CostDescriber::getVertexWeights ( )
inlineprivate
int Isorropia::Epetra::CostDescriber::getGraphEdgeWeights ( int  vertex_global_id,
std::map< int, float > &  wgtMap 
) const
private

getGraphEdgeWeights is called to obtain the graph edge weights for a given vertex (row).

Parameters
vertex_global_idthe global ID of the vertex the caller wants edge weights for
wgtMapa map from the global ID of each vertex neighbor to the weight of the edge formed by the vertex and this neighbor
Returns
the count of the neighbors of vertex_global_id
bool Isorropia::Epetra::CostDescriber::haveGlobalVertexWeights ( ) const
private

haveGlobalVertexWeights returns true if any process in the application has supplied vertex weights, it returns false otherwise.

void Isorropia::Epetra::CostDescriber::setNumGlobalVertexWeights ( int  num)
private

setNumGlobalVertexWeights may be used to set the count of the global number of vertex weights supplied to the CostDescriber

bool Isorropia::Epetra::CostDescriber::haveGlobalGraphEdgeWeights ( ) const
private

haveGlobalGraphEdgeWeights returns true if any process in the application has supplied graph edge weights, it returns false otherwise.

void Isorropia::Epetra::CostDescriber::setNumGlobalGraphEdgeWeights ( int  num)
private

setNumGlobalGraphEdgeWeights may be used to set the count of the global number of graph edge weights supplied to the CostDescriber

bool Isorropia::Epetra::CostDescriber::haveGlobalHypergraphEdgeWeights ( ) const
private

haveGlobalHypergraphEdgeWeights returns true if any process in the application has supplied hyperedge weights, it returns false otherwise.

void Isorropia::Epetra::CostDescriber::setNumGlobalHypergraphEdgeWeights ( int  num)
private

setNumGlobalHypergraphEdgeWeights may be used to set the count of the global number of hyperedge weights supplied to the CostDescriber

void Isorropia::Epetra::CostDescriber::allocate_hg_edge_weights_ ( int  n)
private

Dynamically allocate storage for hypergraph edge weights.

void Isorropia::Epetra::CostDescriber::free_hg_edge_weights_ ( )
private

Free storage used by hypergraph edge weights.

int Isorropia::Epetra::CostDescriber::getEdges ( int  vertexGID,
int  len,
int *  nborGID,
float *  weights 
) const
private

getEdges creates an array of the neighbors and edge weights for given vertex.

Self edges are not included.

Parameters
vertexGIDthe global ID of the vertex (must be one owned by calling process)
lenof preallocated nborGID and weights arrays
nborGIDon return contains the global ID of each vertex neighboring vertexGID, allocated by caller
weightson return contains the weight for each edge formed by the vertices in nborGID
Returns
the number of neighbors in nborGID is returned

Friends And Related Function Documentation

friend class Isorropia::Operator
friend
friend class Isorropia::Epetra::ZoltanLibClass
friend
std::ostream& operator<< ( std::ostream &  ,
const Isorropia::Epetra::CostDescriber cd 
)
friend

Overloaded << operator for CostDescriber object.

Member Data Documentation

Teuchos::RCP<const Epetra_Vector> Isorropia::Epetra::CostDescriber::vertex_weights_
private
Teuchos::RCP<const Epetra_CrsMatrix> Isorropia::Epetra::CostDescriber::graph_edge_weights_
private
std::set<int> Isorropia::Epetra::CostDescriber::graph_self_edges_
private
Teuchos::ParameterList Isorropia::Epetra::CostDescriber::paramlist_
private
int* Isorropia::Epetra::CostDescriber::hg_edge_gids_
private
float* Isorropia::Epetra::CostDescriber::hg_edge_weights_
private
int Isorropia::Epetra::CostDescriber::num_hg_edge_weights_
private
int Isorropia::Epetra::CostDescriber::numGlobalVertexWeights_
private
int Isorropia::Epetra::CostDescriber::numGlobalGraphEdgeWeights_
private
int Isorropia::Epetra::CostDescriber::numGlobalHypergraphEdgeWeights_
private

The documentation for this class was generated from the following file: