Isorropia: Partitioning, Load Balancing and more
Public Member Functions | Private Member Functions | List of all members
Isorropia::CostDescriber Class Referenceabstract

Interface (abstract base class) for describing the weights or costs associated with the vertices and/or edges or hyperedges of the object to be partitioned, ordered or colored. More...

#include <Isorropia_CostDescriber.hpp>

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

Public Member Functions

virtual ~CostDescriber ()
 Destructor. More...
 

Private Member Functions

virtual void setParameters (const Teuchos::ParameterList &paramlist)=0
 Set parameters for the CostDescriber instance. More...
 
virtual bool haveVertexWeights () const =0
 Query whether vertex weights have been supplied by the application. More...
 
virtual int getNumVertices () const =0
 Get the number of vertices for which this process supplied vertex weights. More...
 
virtual void getVertexWeights (int numVertices, int *global_ids, float *weights) const =0
 Get lists of the vertex ids and weights supplied by this process. More...
 
virtual bool haveGraphEdgeWeights () const =0
 Query whether graph edge weights have been supplied by the application. More...
 
virtual int getNumGraphEdges (int vertex_global_id) const =0
 Get the number of graph edges for a specified vertex. More...
 
virtual void getGraphEdgeWeights (int vertex_global_id, int num_neighbors, int *neighbor_global_ids, float *weights) const =0
 Get the graph edge weights for a specified vertex. More...
 
virtual bool haveHypergraphEdgeWeights () const =0
 Query whether hypergraph edge weights have been supplied by the application. More...
 
virtual int getNumHypergraphEdgeWeights () const =0
 Get the number of Hypergraph edges. More...
 
virtual void getHypergraphEdgeWeights (int numEdges, int *global_ids, float *weights) const =0
 Get the hypergraph edge weights that were supplied by this process. More...
 

Detailed Description

Interface (abstract base class) for describing the weights or costs associated with the vertices and/or edges or hyperedges of the object to be partitioned, ordered or colored.

A CostDescriber object is created by the application. If no CostDescriber is supplied by the application, sensible default weights should be used.

Constructor & Destructor Documentation

virtual Isorropia::CostDescriber::~CostDescriber ( )
inlinevirtual

Destructor.

Reimplemented in Isorropia::Epetra::CostDescriber.

Member Function Documentation

virtual void Isorropia::CostDescriber::setParameters ( const Teuchos::ParameterList &  paramlist)
privatepure virtual

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.

Implemented in Isorropia::Epetra::CostDescriber.

virtual bool Isorropia::CostDescriber::haveVertexWeights ( ) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual int Isorropia::CostDescriber::getNumVertices ( ) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual void Isorropia::CostDescriber::getVertexWeights ( int  numVertices,
int *  global_ids,
float *  weights 
) const
privatepure virtual

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.

Implemented in Isorropia::Epetra::CostDescriber.

virtual bool Isorropia::CostDescriber::haveGraphEdgeWeights ( ) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual int Isorropia::CostDescriber::getNumGraphEdges ( int  vertex_global_id) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual void Isorropia::CostDescriber::getGraphEdgeWeights ( int  vertex_global_id,
int  num_neighbors,
int *  neighbor_global_ids,
float *  weights 
) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual bool Isorropia::CostDescriber::haveHypergraphEdgeWeights ( ) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual int Isorropia::CostDescriber::getNumHypergraphEdgeWeights ( ) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.

virtual void Isorropia::CostDescriber::getHypergraphEdgeWeights ( int  numEdges,
int *  global_ids,
float *  weights 
) const
privatepure virtual

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

Implemented in Isorropia::Epetra::CostDescriber.


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