Isorropia: Partitioning, Load Balancing and more
Namespaces | Classes | Functions
Isorropia::Epetra Namespace Reference

Namespaces

 ZoltanLib
 The CostDescriber class describes the vertex, edge and/or hyperedge weights.
 

Classes

class  Colorer
 An implementation of the Colorer interface that operates on Epetra matrices and linear systems. More...
 
class  CostDescriber
 
class  LevelScheduler
 An implementation of the LevelScheduler interface that operates on and Epetra_CrsGraph, representing the non-zeros in a matrix. More...
 
class  Library
 An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...
 
class  Matcher
 An implementation of the Matcher interface that operates on Epetra matrices and Graphs. More...
 
class  Operator
 An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...
 
class  Orderer
 An implementation of the Orderer interface that operates on Epetra matrices and linear systems. More...
 
class  Partitioner
 An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...
 
class  Partitioner2D
 An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...
 
class  Prober
 An implementation of the Prober interface that operates on Epetra matrices and linear systems. More...
 
class  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...
 
class  ZoltanLibClass
 

Functions

Epetra_MultiVector * createBalancedCopy (const Epetra_MultiVector &input_vector)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_MultiVector * createBalancedCopy (const Epetra_MultiVector &input_vector, const Teuchos::ParameterList &paramlist)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_CrsGraph * createBalancedCopy (const Epetra_CrsGraph &input_graph)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_CrsGraph * createBalancedCopy (const Epetra_CrsGraph &input_graph, const Teuchos::ParameterList &paramlist)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_CrsMatrix * createBalancedCopy (const Epetra_CrsMatrix &input_matrix)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_CrsMatrix * createBalancedCopy (const Epetra_CrsMatrix &input_matrix, const Teuchos::ParameterList &paramlist)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_LinearProblem * createBalancedCopy (const Epetra_LinearProblem &input_problem)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
Epetra_LinearProblem * createBalancedCopy (const Epetra_LinearProblem &input_problem, const Teuchos::ParameterList &paramlist)
 createBalancedCopy() creates a copy with a more balanced map. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsMatrix > 
create_balanced_copy (const Epetra_CrsMatrix &input_matrix)
 create_target_map() is an internal function used by Isorropia. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsMatrix > 
create_balanced_copy (const Epetra_CrsMatrix &input_matrix, const Epetra_Vector &row_weights)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsMatrix > 
create_balanced_copy (const Epetra_CrsMatrix &input_matrix, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsMatrix > 
create_balanced_copy (const Epetra_CrsMatrix &input_matrix, CostDescriber &costs, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_RowMatrix > 
create_balanced_copy (const Epetra_RowMatrix &input_matrix)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_RowMatrix > 
create_balanced_copy (const Epetra_RowMatrix &input_matrix, const Epetra_Vector &row_weights)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_RowMatrix > 
create_balanced_copy (const Epetra_RowMatrix &input_matrix, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_RowMatrix > 
create_balanced_copy (const Epetra_RowMatrix &input_matrix, CostDescriber &costs, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsGraph > 
create_balanced_copy (const Epetra_CrsGraph &input_graph)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsGraph > 
create_balanced_copy (const Epetra_CrsGraph &input_matrix, const Epetra_Vector &row_weights)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsGraph > 
create_balanced_copy (const Epetra_CrsGraph &input_graph, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_CrsGraph > 
create_balanced_copy (const Epetra_CrsGraph &input_graph, CostDescriber &costs, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_LinearProblem > 
create_balanced_copy (const Epetra_LinearProblem &input_problem)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_LinearProblem > 
create_balanced_copy (const Epetra_LinearProblem &input_matrix, const Epetra_Vector &row_weights)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_LinearProblem > 
create_balanced_copy (const Epetra_LinearProblem &input_problem, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_LinearProblem > 
create_balanced_copy (const Epetra_LinearProblem &input_problem, CostDescriber &costs, const Teuchos::ParameterList &paramlist)
 create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix. More...
 
__deprecated Teuchos::RCP
< Epetra_MultiVector > 
create_balanced_copy (const Epetra_MultiVector &coords, const Teuchos::ParameterList &paramlist)
 
__deprecated Teuchos::RCP
< Epetra_MultiVector > 
create_balanced_copy (const Epetra_MultiVector &coords, const Epetra_MultiVector &weights, const Teuchos::ParameterList &paramlist)
 
__deprecated Teuchos::RCP
< Epetra_MultiVector > 
create_balanced_copy (const Epetra_MultiVector &coords)
 
__deprecated Teuchos::RCP
< Epetra_MultiVector > 
create_balanced_copy (const Epetra_MultiVector &coords, const Epetra_MultiVector &weights)
 
Teuchos::RCP< Epetra_CrsMatrix > redistribute_rows (const Epetra_CrsMatrix &input_matrix, const Epetra_Map &target_rowmap, Epetra_Import *importer=0)
 redistribute_rows() is an internal Isorropia function, not part of the API. More...
 
Teuchos::RCP< Epetra_CrsMatrix > redistribute_rows (const Epetra_RowMatrix &input_matrix, const Epetra_Map &target_rowmap, Epetra_Import *importer=0)
 redistribute_rows() is an internal Isorropia function, not part of the API. More...
 
Teuchos::RCP< Epetra_CrsGraph > redistribute_rows (const Epetra_CrsGraph &input_graph, const Epetra_Map &target_rowmap, Epetra_Import *importer=0)
 Return a new Epetra_CrsGraph object constructed with target_rowmap, and with the contents of input_graph imported into it. More...
 
Teuchos::RCP< Epetra_MultiVector > redistribute (const Epetra_MultiVector &input, const Epetra_BlockMap &target_map, Epetra_Import *importer=0)
 Return a new Epetra_MultiVector object constructed with target_map, and with the contents of 'input' imported into it. More...
 
Teuchos::RCP< Epetra_Vector > redistribute (const Epetra_Vector &input, const Epetra_Map &target_map, Epetra_Import *importer=0)
 Return a new Epetra_Vector object constructed with target_map, and with the contents of 'input' imported into it. More...
 

Function Documentation

Epetra_MultiVector* Isorropia::Epetra::createBalancedCopy ( const Epetra_MultiVector &  input_vector)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Examples:
matrix_1.cpp.
Epetra_MultiVector* Isorropia::Epetra::createBalancedCopy ( const Epetra_MultiVector &  input_vector,
const Teuchos::ParameterList &  paramlist 
)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_CrsGraph* Isorropia::Epetra::createBalancedCopy ( const Epetra_CrsGraph &  input_graph)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_CrsGraph* Isorropia::Epetra::createBalancedCopy ( const Epetra_CrsGraph &  input_graph,
const Teuchos::ParameterList &  paramlist 
)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_CrsMatrix* Isorropia::Epetra::createBalancedCopy ( const Epetra_CrsMatrix &  input_matrix)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_CrsMatrix* Isorropia::Epetra::createBalancedCopy ( const Epetra_CrsMatrix &  input_matrix,
const Teuchos::ParameterList &  paramlist 
)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_LinearProblem* Isorropia::Epetra::createBalancedCopy ( const Epetra_LinearProblem &  input_problem)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

Epetra_LinearProblem* Isorropia::Epetra::createBalancedCopy ( const Epetra_LinearProblem &  input_problem,
const Teuchos::ParameterList &  paramlist 
)

createBalancedCopy() creates a copy with a more balanced map.

The caller should free the copy after use.

__deprecated Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsMatrix &  input_matrix)

create_target_map() is an internal function used by Isorropia.

Given a Partitioner object, create a target map representing the new partitioning. This function calls partitioner.compute_partitioning() if it has not already been called. create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to balance row weights while minimizing cuts in the hyperedges.

The non-Zoltan rebalancing is 1-D, row-wise, and attempts to make the number of nonzeros equal in each partition. I.e., it is equivalent to specifying a weighted rebalance where the weights are the number of nonzeros in each row.

__deprecated Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsMatrix &  input_matrix,
const Epetra_Vector &  row_weights 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

The balancing algorithm (Zoltan or non-Zoltan) will use the row weights provided for the matrix rows.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign a unit weight to each column (hyperedge). It will attempt to balance row weights while minimizing cuts in the hyperedges.

The non-Zoltan rebalancing is a quick 1-D, row-wise balancing.

__deprecated Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsMatrix &  input_matrix,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

__deprecated Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsMatrix &  input_matrix,
CostDescriber &  costs,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has. It ignores the costs object.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan. The weights provided in the CostDescriber object will be supplied to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method.

If an empty CostDescriber is supplied, we will assign unit weight to each row (vertex), and unit weight to each column (hyperedge), in the case of hypergraph partitioning. We will assign unit weight to each row (vertex) and each non zero (edge) in the case of graph partitioning.

__deprecated Teuchos::RCP<Epetra_RowMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_RowMatrix &  input_matrix)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

The non-Zoltan rebalancing is 1-D, row-wise, and attempts to make the number of nonzeros equal in each partition. I.e., it is equivalent to specifying a weighted rebalance where the weights are the number of nonzeros in each row.

__deprecated Teuchos::RCP<Epetra_RowMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_RowMatrix &  input_matrix,
const Epetra_Vector &  row_weights 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

The balancing algorithm (Zoltan or non-Zoltan) will use the row weights provided for the matrix rows.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign a unit weight to each column (hyperedge). It will attempt to balance row weights while minimizing cuts in the hyperedges.

The non-Zoltan rebalancing is a quick 1-D, row-wise balancing.

__deprecated Teuchos::RCP<Epetra_RowMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_RowMatrix &  input_matrix,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

__deprecated Teuchos::RCP<Epetra_RowMatrix> Isorropia::Epetra::create_balanced_copy ( const Epetra_RowMatrix &  input_matrix,
CostDescriber &  costs,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has. The costs object is ignored in this case.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan. The weights provided in the CostDescriber object will be supplied to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method.

If an empty CostDescriber is supplied, we will assign unit weight to each row (vertex), and unit weight to each column (hyperedge), in the case of hypergraph partitioning. We will assign unit weight to each row (vertex) and each non zero (edge) in the case of graph partitioning.

__deprecated Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsGraph &  input_graph)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

The non-Zoltan rebalancing is 1-D, row-wise, and attempts to make the number of nonzeros equal in each partition. I.e., it is equivalent to specifying a weighted rebalance where the weights are the number of nonzeros in each row.

__deprecated Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsGraph &  input_matrix,
const Epetra_Vector &  row_weights 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

The balancing algorithm (Zoltan or non-Zoltan) will use the row weights provided for the matrix rows.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign a unit weight to each column (hyperedge). It will attempt to balance row weights while minimizing cuts in the hyperedges.

The non-Zoltan rebalancing is a quick 1-D, row-wise balancing.

__deprecated Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsGraph &  input_graph,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

__deprecated Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::create_balanced_copy ( const Epetra_CrsGraph &  input_graph,
CostDescriber &  costs,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan. The weights provided in the CostDescriber object will be supplied to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method.

If an empty CostDescriber is supplied, we will assign unit weight to each row (vertex), and unit weight to each column (hyperedge), in the case of hypergraph partitioning. We will assign unit weight to each row (vertex) and each non zero (edge) in the case of graph partitioning.

__deprecated Teuchos::RCP<Epetra_LinearProblem> Isorropia::Epetra::create_balanced_copy ( const Epetra_LinearProblem &  input_problem)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

The non-Zoltan rebalancing is 1-D, row-wise, and attempts to make the number of nonzeros equal in each partition. I.e., it is equivalent to specifying a weighted rebalance where the weights are the number of nonzeros in each row.

__deprecated Teuchos::RCP<Epetra_LinearProblem> Isorropia::Epetra::create_balanced_copy ( const Epetra_LinearProblem &  input_matrix,
const Epetra_Vector &  row_weights 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia will use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

The balancing algorithm (Zoltan or non-Zoltan) will use the row weights provided for the matrix rows.

Zoltan will perform hypergraph partitioning using its own PHG method. It will assign a unit weight to each column (hyperedge). It will attempt to balance row weights while minimizing cuts in the hyperedges.

The non-Zoltan rebalancing is a quick 1-D, row-wise balancing.

__deprecated Teuchos::RCP<Epetra_LinearProblem> Isorropia::Epetra::create_balanced_copy ( const Epetra_LinearProblem &  input_problem,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method. It will assign unit weight to each row (vertex), and unit weight to each column (hyperedge). It will attempt to minimize cuts in the hyperedges.

__deprecated Teuchos::RCP<Epetra_LinearProblem> Isorropia::Epetra::create_balanced_copy ( const Epetra_LinearProblem &  input_problem,
CostDescriber &  costs,
const Teuchos::ParameterList &  paramlist 
)

create_balanced_copy(), which is part of the Isorropia API, is used to create and return a balanced copy of an input Epetra_CrsMatrix.

Isorropia can use Zoltan to perform partitioning if it has been configured with Zoltan support. If Zoltan is not available, it will perform its own simple row partitioning.

In addition, if the parameter "PARTITIONING_METHOD" is set to "SIMPLE_LINEAR", Isorropia will perform its own simple row partitioning, even if Zoltan is available. This method balances the rows across processes after assigning each row a weight equal to the number of non zeros it has. This costs object is ignored in this case.

If Zoltan is called to do the partitioning, any parameters set in a "Zoltan" sublist of the paramlist are provided to Zoltan. The weights provided in the CostDescriber object will be supplied to Zoltan.

Refer to the Zoltan users guide for specific parameters that Zoltan recognizes. A couple of important ones are "LB_METHOD" (valid values include "GRAPH", "HYPERGRAPH"), "DEBUG_LEVEL" (valid values are 0 to 10, default is 1), etc.

If no Zoltan parameters are provided, Zoltan will perform hypergraph partitioning using its own PHG method.

If an empty CostDescriber is supplied, we will assign unit weight to each row (vertex), and unit weight to each column (hyperedge), in the case of hypergraph partitioning. We will assign unit weight to each row (vertex) and each non zero (edge) in the case of graph partitioning.

__deprecated Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::create_balanced_copy ( const Epetra_MultiVector &  coords,
const Teuchos::ParameterList &  paramlist 
)
__deprecated Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::create_balanced_copy ( const Epetra_MultiVector &  coords,
const Epetra_MultiVector &  weights,
const Teuchos::ParameterList &  paramlist 
)
__deprecated Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::create_balanced_copy ( const Epetra_MultiVector &  coords)
__deprecated Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::create_balanced_copy ( const Epetra_MultiVector &  coords,
const Epetra_MultiVector &  weights 
)
Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::redistribute_rows ( const Epetra_CrsMatrix &  input_matrix,
const Epetra_Map &  target_rowmap,
Epetra_Import *  importer = 0 
)

redistribute_rows() is an internal Isorropia function, not part of the API.

Return a new Epetra_CrsMatrix object constructed with target_rowmap, and with the contents of input_matrix imported into it.

The caller is responsible for deleting the returned object.

param input_matrix Source/input object.

param target_rowmap Target rowmap, required to be compatible with input_matrix.RowMap() in terms of number-of-elements, etc.

param importer Optional argument. If importer is supplied, it will be used to perform the import operation. Otherwise, a temporary importer will be created and used.

Teuchos::RCP<Epetra_CrsMatrix> Isorropia::Epetra::redistribute_rows ( const Epetra_RowMatrix &  input_matrix,
const Epetra_Map &  target_rowmap,
Epetra_Import *  importer = 0 
)

redistribute_rows() is an internal Isorropia function, not part of the API.

Return a new Epetra_CrsMatrix object constructed with target_rowmap, and with the contents of input_matrix imported into it.

The caller is responsible for deleting the returned object.

param input_matrix Source/input object.

param target_rowmap Target rowmap, required to be compatible with input_matrix.RowMatrixRowMap() in terms of number-of-elements, etc.

param importer Optional argument. If importer is supplied, it will be used to perform the import operation. Otherwise, a temporary importer will be created and used.

Teuchos::RCP<Epetra_CrsGraph> Isorropia::Epetra::redistribute_rows ( const Epetra_CrsGraph &  input_graph,
const Epetra_Map &  target_rowmap,
Epetra_Import *  importer = 0 
)

Return a new Epetra_CrsGraph object constructed with target_rowmap, and with the contents of input_graph imported into it.

param input_graph Source/input object.

param target_rowmap Target rowmap, required to be compatible with input_graph.RowMap() in terms of number-of-elements, etc.

param importer Optional argument. If importer is supplied, it will be used to perform the import operation. Otherwise, a temporary importer will be created and used.

Teuchos::RCP<Epetra_MultiVector> Isorropia::Epetra::redistribute ( const Epetra_MultiVector &  input,
const Epetra_BlockMap &  target_map,
Epetra_Import *  importer = 0 
)

Return a new Epetra_MultiVector object constructed with target_map, and with the contents of 'input' imported into it.

param input Source/input object.

param target_map Target map, required to be compatible with input.Map() in terms of number-of-elements, etc.

param importer Optional argument. If importer is supplied, it will be used to perform the import operation. Otherwise, a temporary importer will be created and used.

Teuchos::RCP<Epetra_Vector> Isorropia::Epetra::redistribute ( const Epetra_Vector &  input,
const Epetra_Map &  target_map,
Epetra_Import *  importer = 0 
)

Return a new Epetra_Vector object constructed with target_map, and with the contents of 'input' imported into it.

param input Source/input object.

param target_map Target map, required to be compatible with input.RowMap() in terms of number-of-elements, etc.

param importer Optional argument. If importer is supplied, it will be used to perform the import operation. Otherwise, a temporary importer will be created and used.