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

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

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

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 (const Epetra_CrsGraph *inputGraph, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_CrsGraph *inputGraph, CostDescriber *costs, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_RowMatrix *inputMatrix, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_RowMatrix *inputMatrix, CostDescriber *costs, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_MultiVector *coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_MultiVector *coords, const Epetra_MultiVector *weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_BlockMap *inputMap, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_CrsGraph *inputGraph, const Epetra_MultiVector *coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_CrsGraph *inputGraph, CostDescriber *costs, const Epetra_MultiVector *coords, const Epetra_MultiVector *weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_RowMatrix *inputMatrix, const Epetra_MultiVector *coords, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
 Isorropia::Epetra::Partitioner::Partitioner (const Epetra_RowMatrix *inputMatrix, CostDescriber *costs, const Epetra_MultiVector *coords, const Epetra_MultiVector *weights, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_partitioning_now=true)
 
void Isorropia::Epetra::Partitioner::createNewMap (Epetra_Map *&outputMap)
 Create a new Epetra_Map corresponding to the new partition. More...
 
 Isorropia::Epetra::Redistributor::Redistributor (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 (Epetra_Map *target_map)
 This constructor sets the target map for the redistribution. More...
 
void Isorropia::Epetra::Redistributor::redistribute (const Epetra_CrsGraph &input_graph, Epetra_CrsGraph *&outputGraphPtr, bool callFillComplete=true)
 Method to accept a Epetra_CrsGraph object, and return a redistributed Epetra_CrsGraph object. More...
 
void Isorropia::Epetra::Redistributor::redistribute (const Epetra_CrsMatrix &inputMatrix, Epetra_CrsMatrix *&outputMatrix, bool callFillComplete=true)
 Method to accept a Epetra_CrsMatrix object, and return a redistributed Epetra_CrsMatrix object. More...
 
void Isorropia::Epetra::Redistributor::redistribute (const Epetra_RowMatrix &inputMatrix, Epetra_CrsMatrix *&outputMatrix, bool callFillComplete=true)
 Method to accept a Epetra_RowMatrix object, and return a redistributed Epetra_CrsMatrix object. More...
 
void Isorropia::Epetra::Redistributor::redistribute (const Epetra_Vector &inputVector, Epetra_Vector *&outputVector)
 Method to accept a Epetra_Vector object, and return a redistributed Epetra_Vector object. More...
 
void Isorropia::Epetra::Redistributor::redistribute (const Epetra_MultiVector &inputVector, Epetra_MultiVector *&outputVector)
 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 Pointers.

Partitioning and Load-balancing with Pointers

One interface to the partitioning in Isorropia is to pass a pointer to the Epetra object to the Isorropia partitioner. This interface is not as safe as the RCP interface, and you need to be careful not to delete the Epetra object when the pointer is still being used inside the partitioner. The following example demonstrates how to partition a matrix with hypergraph partitioning using the pointer interface.

  // loadMatrix() allocates and reads in an Epetra_CrsMatrix from a file
  const Epetra_CrsMatrix *A = 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(A,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 pointer interface to the redistributor accepts a pointer to an Isorropia::Epetra::Partitioner object and and sets a pointer to the newly allocated matrix/graph/etc. Below we demonstrate how the redistributor works with the pointer interface.

  // loadMatrix() allocates and reads in an Epetra_RowMatrix from a file
  Epetra_CrsMatrix *A; = loadMatrix(filename);
  Epetra_CrsMatrix *newMatrix;

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

  Isorropia::Epetra::Partitioner *partitioner = new Isorropia::Epetra::Partitioner(A, params,false);
  partitioner->partition();

  Isorropia::Epetra::Redistributor rd(partitioner);
  rd.redistribute(*rowmatrix,newMatrix,true);

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

Function Documentation

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

Create a new Epetra_Map corresponding to the new partition.

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

Parameters
[out]outputMapEpetra_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 ( 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 ( 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
void Isorropia::Epetra::Redistributor::redistribute ( const Epetra_CrsGraph &  input_graph,
Epetra_CrsGraph *&  outputGraphPtr,
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.
outputGraphPtr(out) pointer to the new redistributed graph
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.
void Isorropia::Epetra::Redistributor::redistribute ( const Epetra_CrsMatrix &  inputMatrix,
Epetra_CrsMatrix *&  outputMatrix,
bool  callFillComplete = true 
)

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

Parameters
inputMatrix(in) the matrix for which we want a new matrix that is distributed according to the partitioner with which this Redistributor was created.
outputMatrix(out) pointer to the new redistributed matrix
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.
void Isorropia::Epetra::Redistributor::redistribute ( const Epetra_RowMatrix &  inputMatrix,
Epetra_CrsMatrix *&  outputMatrix,
bool  callFillComplete = true 
)

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

Parameters
inputMatrix(in) the row matrix for which we want a new matrix that is distributed according to the partitioner with which this Redistributor was created.
outputMatrix(out) pointer to the new redistributed matrix
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.
void Isorropia::Epetra::Redistributor::redistribute ( const Epetra_Vector &  inputVector,
Epetra_Vector *&  outputVector 
)

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

Parameters
inputVector(in) the vector for which we want a new vector that is distributed according to the partitioner with which this Redistributor was created.
outputVector(out) pointer to the new redistributed vector
void Isorropia::Epetra::Redistributor::redistribute ( const Epetra_MultiVector &  inputVector,
Epetra_MultiVector *&  outputVector 
)

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

Parameters
inputVector(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.
outputVector(out) a reference counted pointer to the new redistributed multi vector