Isorropia: Partitioning, Load Balancing and more
Classes | Functions
Isorropia Partitioning and Load Balancing Methods with RCPs

Here we describe the partitioning interface that uses RCPs. More...

Collaboration diagram for Isorropia Partitioning and Load Balancing Methods with RCPs:

Classes

class  Isorropia::Epetra::Partitioner
 An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...
 
class  Isorropia::Epetra::Redistributor
 Class which is constructed with a Partitioner instance, and provides several methods for redistributing Epetra objects given the partitioning computed by the Partitioner object. More...
 

Functions

 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_CrsGraph > inputGraph, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_CrsGraph > inputGraph, Teuchos::RCP< CostDescriber > costs, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_RowMatrix > inputMatrix, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_RowMatrix > inputMatrix, Teuchos::RCP< CostDescriber > costs, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_MultiVector > coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_BlockMap > inputMap, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_CrsGraph > inputGraph, Teuchos::RCP< const Epetra_MultiVector > coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_CrsGraph > inputGraph, Teuchos::RCP< CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_RowMatrix > inputMatrix, Teuchos::RCP< const Epetra_MultiVector > coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (Teuchos::RCP< const Epetra_RowMatrix > inputMatrix, Teuchos::RCP< CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
Teuchos::RCP< Epetra_Map > Isorropia::Epetra::Partitioner::createNewMap ()
 Create a new Epetra_Map corresponding to the new partition. More...
 
 Isorropia::Epetra::Redistributor::Redistributor (Teuchos::RCP< Isorropia::Epetra::Partitioner > partitioner)
 This constructor calls the Isorropia::Epetra::Partitioner::partition method on the partitioner if it has not already been called. More...
 
 Isorropia::Epetra::Redistributor::Redistributor (Teuchos::RCP< Epetra_Map > target_map)
 This constructor sets the target map for the redistribution. More...
 
Teuchos::RCP< Epetra_CrsGraph > Isorropia::Epetra::Redistributor::redistribute (const Epetra_CrsGraph &input_graph, bool callFillComplete=true)
 Method to accept a Epetra_CrsGraph object, and return a redistributed Epetra_CrsGraph object. More...
 
Teuchos::RCP< Epetra_CrsMatrix > Isorropia::Epetra::Redistributor::redistribute (const Epetra_CrsMatrix &input_matrix, bool callFillComplete=true)
 Method to accept a Epetra_CrsMatrix object, and return a redistributed Epetra_CrsMatrix object. More...
 
Teuchos::RCP< Epetra_CrsMatrix > Isorropia::Epetra::Redistributor::redistribute (const Epetra_RowMatrix &input_matrix, bool callFillComplete=true)
 Method to accept a Epetra_RowMatrix object, and return a redistributed Epetra_CrsMatrix object. More...
 
Teuchos::RCP< Epetra_Vector > Isorropia::Epetra::Redistributor::redistribute (const Epetra_Vector &input_vector)
 Method to accept a Epetra_Vector object, and return a redistributed Epetra_Vector object. More...
 
Teuchos::RCP< Epetra_MultiVector > Isorropia::Epetra::Redistributor::redistribute (const Epetra_MultiVector &input_vector)
 Method to accept a Epetra_MultiVector object, and return a redistributed Epetra_MultiVector object. More...
 

Detailed Description

Here we describe the partitioning interface that uses RCPs.

Partitioning and Load-balancing with RCPs

The safer interface to the partitioning in Isorropia is to use the Teuchos reference counter pointer (RCP) and pass a RCP to the Epetra object to the Isorropia partitioner. We suggest using interface in order to prevent memory leaks if you are familiar with RCPs and willing to learn the basics of Teuchos::RCP. The following example demonstrates how to partition a matrix with hypergraph partitioning using the RCP interface.

  // loadMatrix() allocates and reads in an Epetra_RowMatrix from a file
  Teuchos::RCP<const Epetra_RowMatrix> rowmatrix = loadMatrix(filename);

  Teuchos::ParameterList params;
  params.set("PARTITIONING_METHOD", "HYPERGRAPH");  
  params.set("BALANCE OBJECTIVE","NONZEROS");
  params.set("IMBALANCE TOL","1.03");

  Isorropia::Epetra::Partitioner partitioner(rowmatrix,params,false);
  partitioner.partition();

After partitioning the matrix, the Isorropia::Epetra::Redistributor may be used to return a new matrix that contains the data of the original matrix/graph/etc. but redistributed based on a partition. The safer RCP interface to the redistributor accepts an RCP to an Isorropia::Epetra::Partitioner object and returns an RCP to a new matrix/graph/etc. Below we demonstrate how the redistributor works with the RCP interface.

  // loadMatrix() allocates and reads in an Epetra_RowMatrix from a file
  Teuchos::RCP<const Epetra_RowMatrix> rowmatrix = loadMatrix(filename);

  Teuchos::ParameterList params;
  params.set("IMBALANCE TOL","1.03");
  params.set("BALANCE OBJECTIVE","NONZEROS");
  params.set("PARTITIONING_METHOD", "HYPERGRAPH");

  Teuchos::RCP<Isorropia::Epetra::Partitioner> partitioner =
       Teuchos::rcp(new Isorropia::Epetra::Partitioner(rowmatrix, params,false));
  partitioner->partition();

  Isorropia::Epetra::Redistributor rd(partitioner);
  Teuchos::RCP<Epetra_CrsMatrix> newmatrix = rd.redistribute(*rowmatrix,true);

Additional information and functionality for partitioning, in general, is documented here:

Function Documentation

Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_CrsGraph >  inputGraph,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_CrsGraph >  inputGraph,
Teuchos::RCP< CostDescriber costs,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_RowMatrix >  inputMatrix,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_RowMatrix >  inputMatrix,
Teuchos::RCP< CostDescriber costs,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_MultiVector >  coords,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_BlockMap >  inputMap,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_CrsGraph >  inputGraph,
Teuchos::RCP< const Epetra_MultiVector >  coords,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_CrsGraph >  inputGraph,
Teuchos::RCP< CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_RowMatrix >  inputMatrix,
Teuchos::RCP< const Epetra_MultiVector >  coords,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Isorropia::Epetra::Partitioner::Partitioner ( Teuchos::RCP< const Epetra_RowMatrix >  inputMatrix,
Teuchos::RCP< CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_partitioning_now = true 
)
Teuchos::RCP<Epetra_Map> Isorropia::Epetra::Partitioner::createNewMap ( )

Create a new Epetra_Map corresponding to the new partition.

This method is essentially used by the Isorropia::Epetra::Redistributor object.

Returns
Epetra_Map that contains the new distribution of elements.
Precondition
The number of parts might be the same or lower than the number of processors.
Isorropia::Epetra::Redistributor::Redistributor ( Teuchos::RCP< Isorropia::Epetra::Partitioner partitioner)

This constructor calls the Isorropia::Epetra::Partitioner::partition method on the partitioner if it has not already been called.

Parameters
partitioner(in) this input partitioner determines the new partitioning to be created when Isorropia::Epetra::Redistributor::redistribute is called
Isorropia::Epetra::Redistributor::Redistributor ( Teuchos::RCP< Epetra_Map >  target_map)

This constructor sets the target map for the redistribution.

Parameters
target_map(in) this input map determines the new matrices/vectors to be created when Isorropia::Epetra::Redistributor::redistribute is called
Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::Redistributor::redistribute ( const Epetra_CrsGraph &  input_graph,
bool  callFillComplete = true 
)

Method to accept a Epetra_CrsGraph object, and return a redistributed Epetra_CrsGraph object.

Parameters
input_graph(in) the graph for which we want a new graph that is distributed according to the partitioner with which this Redistributor was created.
callFillComplete(in) The new graph is FillComplete'd if callFillComplete is true. In that case, the range map is set to equal the row map. The domain map will equal the range map, unless the input_graph has different domain and range maps, in which case the original domain map is preserved. By default callFillComplete is true.
Returns
a reference counted pointer to the new redistributed graph
Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::Redistributor::redistribute ( const Epetra_CrsMatrix &  input_matrix,
bool  callFillComplete = true 
)

Method to accept a Epetra_CrsMatrix object, and return a redistributed Epetra_CrsMatrix object.

Parameters
input_matrix(in) the matrix for which we want a new matrix that is distributed according to the partitioner with which this Redistributor was created.
callFillComplete(in) The new matrix is FillComplete'd if callFillComplete is true. In that case, the range map is set to equal the row map. The domain map will equal the range map, unless the input_matrix has different domain and range maps, in which case the original domain map is preserved. By default callFillComplete is true.
Returns
a reference counted pointer to the new redistributed matrix
Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::Redistributor::redistribute ( const Epetra_RowMatrix &  input_matrix,
bool  callFillComplete = true 
)

Method to accept a Epetra_RowMatrix object, and return a redistributed Epetra_CrsMatrix object.

Parameters
input_matrix(in) the row matrix for which we want a new matrix that is distributed according to the partitioner with which this Redistributor was created.
callFillComplete(in) The new matrix is FillComplete'd if callFillComplete is true. In that case, the range map is set to equal the row map. The domain map will equal the range map, unless the input_matrix has different domain and range maps, in which case the original domain map is preserved. By default callFillComplete is true.
Returns
a reference counted pointer to the new redistributed matrix
Teuchos::RCP<Epetra_Vector> Isorropia::Epetra::Redistributor::redistribute ( const Epetra_Vector &  input_vector)

Method to accept a Epetra_Vector object, and return a redistributed Epetra_Vector object.

Parameters
input_vector(in) the vector for which we want a new vector that is distributed according to the partitioner with which this Redistributor was created.
Returns
a reference counted pointer to the new redistributed vector
Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::Redistributor::redistribute ( const Epetra_MultiVector &  input_vector)

Method to accept a Epetra_MultiVector object, and return a redistributed Epetra_MultiVector object.

Parameters
input_vector(in) the multi vector for which we want a new multi vector that is distributed according to the partitioner with which this Redistributor was created.
Returns
a reference counted pointer to the new redistributed multi vector