Isorropia: Partitioning, Load Balancing and more
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
Isorropia::Epetra::Operator Class Referenceabstract

An implementation of the Partitioner interface that operates on Epetra matrices and linear systems. More...

#include <Isorropia_EpetraOperator.hpp>

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

Public Member Functions

 Operator (Teuchos::RCP< const Epetra_CrsGraph > input_graph, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_CrsGraph object. More...
 
 Operator (Teuchos::RCP< const Epetra_BlockMap > input_map, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_BlockMap object. More...
 
 Operator (Teuchos::RCP< const Epetra_CrsGraph > input_graph, Teuchos::RCP< const Epetra_MultiVector > input_coords, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_CrsGraph object and an Epetra_MultiVector object of coordinates. More...
 
 Operator (Teuchos::RCP< const Epetra_CrsGraph > input_graph, Teuchos::RCP< CostDescriber > costs, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_CrsGraph object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights. More...
 
 Operator (Teuchos::RCP< const Epetra_CrsGraph > input_graph, Teuchos::RCP< CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_CrsGraph object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights. More...
 
 Operator (Teuchos::RCP< const Epetra_RowMatrix > input_matrix, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_RowMatrix object. More...
 
 Operator (Teuchos::RCP< const Epetra_RowMatrix > input_matrix, Teuchos::RCP< const Epetra_MultiVector > coords, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_RowMatrix object and an Epetra_MultiVector object of coordinates. More...
 
 Operator (Teuchos::RCP< const Epetra_RowMatrix > input_matrix, Teuchos::RCP< CostDescriber > costs, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_RowMatrix object and a CostDescriber. More...
 
 Operator (Teuchos::RCP< const Epetra_RowMatrix > input_matrix, Teuchos::RCP< CostDescriber > costs, Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_RowMatrix object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights. More...
 
 Operator (Teuchos::RCP< const Epetra_MultiVector > coords, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_MultiVector object and a ParameterList. More...
 
 Operator (Teuchos::RCP< const Epetra_MultiVector > coords, Teuchos::RCP< const Epetra_MultiVector > weights, const Teuchos::ParameterList &paramlist, int base)
 Constructor that accepts an Epetra_MultiVector object and a ParameterList. More...
 
virtual ~Operator ()
 Destructor. More...
 
void setParameters (const Teuchos::ParameterList &paramlist)
 setParameters() is an internal method which handles the parameters from a Teuchos::ParameterList object. More...
 
Teuchos::RCP
< Isorropia::Epetra::CostDescriber > & 
getCosts ()
 Get the cost object. More...
 
virtual void compute (bool force_compute)=0
 Method which does the work of computing a new partitioning/coloring/ordering, depending on the child class used. More...
 
bool alreadyComputed () const
 Query whether compute_operation() has already been called. More...
 
int numProperties () const
 Return the number of different values used for "properties". More...
 
int numLocalProperties () const
 Return the number of different values used for "properties" for this process only. More...
 
virtual const int & operator[] (int myElem) const
 Return the new partition ID for a given element that resided locally in the old operation. More...
 
virtual int numElemsWithProperty (int property) const
 Return the number of elements in a given partition. More...
 
virtual void elemsWithProperty (int property, int *elementList, int len) const
 Fill user-allocated list (of length len) with the global element ids to be located in the given partition. More...
 
virtual int extractPropertiesCopy (int len, int &size, int *array) const
 Copy a part of the property array. More...
 
virtual int extractPropertiesView (int &size, const int *&array) const
 Give access of the property array that is owned by the current processor. More...
 

Protected Member Functions

void computeNumberOfProperties ()
 

Protected Attributes

Teuchos::RCP< const
Epetra_BlockMap > 
input_map_
 
Teuchos::RCP< const
Epetra_CrsGraph > 
input_graph_
 
Teuchos::RCP< const
Epetra_RowMatrix > 
input_matrix_
 
Teuchos::RCP< const
Epetra_MultiVector > 
input_coords_
 
Teuchos::RCP
< Isorropia::Epetra::CostDescriber
costs_
 
Teuchos::RCP< const
Epetra_MultiVector > 
weights_
 
Teuchos::ParameterList paramlist_
 
int exportsSize_
 
std::vector< int > imports_
 
std::vector< int > properties_
 
bool operation_already_computed_
 
int global_num_vertex_weights_
 
int global_num_graph_edge_weights_
 
int global_num_hg_edge_weights_
 
Teuchos::RCP< Librarylib_
 
int base_
 

Private Member Functions

void paramsToUpper (Teuchos::ParameterList &, int &changed, bool rmUnderscore=true)
 
void stringToUpper (std::string &s, int &changed, bool rmUnderscore=false)
 

Private Attributes

int numberOfProperties_
 
int localNumberOfProperties_
 
std::vector< int > numberElemsByProperties_
 

Detailed Description

An implementation of the Partitioner interface that operates on Epetra matrices and linear systems.

Constructor & Destructor Documentation

Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_CrsGraph >  input_graph,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_CrsGraph object.

Parameters
input_graphMatrix-graph object for which a new operation is to be computed.
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseindex base

If the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.

Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_BlockMap >  input_map,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_BlockMap object.

Parameters
[in]input_mapBlockMap object for which a new operation is to be computed.
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_CrsGraph >  input_graph,
Teuchos::RCP< const Epetra_MultiVector >  input_coords,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_CrsGraph object and an Epetra_MultiVector object of coordinates.

Parameters
input_graphMatrix-graph object for which a new operation is to be computed.
[in]input_coordsThe input geometric coordinates (represented in a multivector)
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_CrsGraph >  input_graph,
Teuchos::RCP< CostDescriber costs,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_CrsGraph object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights.

Parameters
input_graphMatrix-graph object for which a new operation is to be computed.
costsCostDescriber object which allows for user-specified weights
paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_CrsGraph >  input_graph,
Teuchos::RCP< CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_CrsGraph object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights.

Parameters
input_graphMatrix-graph object for which a new operation is to be computed.
costsCostDescriber object which allows for user-specified weights
coordsThe input geometric coordinates (represented in a multivector)
weightsA one or more dimensional weight for each of the The input geometric coordinates
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_RowMatrix >  input_matrix,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_RowMatrix object.

Parameters
input_matrixMatrix object for which a new operation is to be computed.
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIndex base

If the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.

Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_RowMatrix >  input_matrix,
Teuchos::RCP< const Epetra_MultiVector >  coords,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_RowMatrix object and an Epetra_MultiVector object of coordinates.

Parameters
input_matrixMatrix object for which a new operation is to be computed.
coordsThe input geometric coordinates (represented in a multivector)
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_RowMatrix >  input_matrix,
Teuchos::RCP< CostDescriber costs,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_RowMatrix object and a CostDescriber.

Parameters
input_matrixMatrix object for which a new operation is to be computed.
costsCostDescriber object which allows for user-specified weights
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_RowMatrix >  input_matrix,
Teuchos::RCP< CostDescriber costs,
Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_RowMatrix object, a CostDescriber, an Epetra_MultiVector object of coordinates, and an Epetra_MultiVector object of coordinate weights.

Parameters
input_matrixMatrix object for which a new operation is to be computed.
costsCostDescriber object which allows for user-specified weights
coordsThe input geometric coordinates (represented in a multivector)
weightsA one or more dimensional weight for each of the The input geometric coordinates
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseIf the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.
Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_MultiVector >  coords,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_MultiVector object and a ParameterList.

Parameters
coordsThe input geometric coordinates (represented in a multivector)
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseindex base

If the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.

Isorropia::Epetra::Operator::Operator ( Teuchos::RCP< const Epetra_MultiVector >  coords,
Teuchos::RCP< const Epetra_MultiVector >  weights,
const Teuchos::ParameterList &  paramlist,
int  base 
)

Constructor that accepts an Epetra_MultiVector object and a ParameterList.

Parameters
coordsThe input geometric coordinates (represented in a multivector)
weightsA one or more dimensional weight for each of the The input geometric coordinates
[in]paramlistTeuchos::ParameterList which will be copied to an internal ParameterList attribute.
[in]baseindex base

If the ParameterList object contains a sublist named "Zoltan", then the Zoltan library is used to perform the operation. Also, any parameters in the "Zoltan" sublist will be relayed directly to Zoltan.

virtual Isorropia::Epetra::Operator::~Operator ( )
virtual

Destructor.

Reimplemented from Isorropia::Operator.

Member Function Documentation

void Isorropia::Epetra::Operator::setParameters ( const Teuchos::ParameterList &  paramlist)
virtual

setParameters() is an internal method which handles the parameters from a Teuchos::ParameterList object.

Implements Isorropia::Operator.

Teuchos::RCP<Isorropia::Epetra::CostDescriber>& Isorropia::Epetra::Operator::getCosts ( )
inline

Get the cost object.

virtual void Isorropia::Epetra::Operator::compute ( bool  forceRecomputing)
pure virtual

Method which does the work of computing a new partitioning/coloring/ordering, depending on the child class used.

Parameters
forceRecomputingOptional argument defaults to false. Depending on the implementation, compute() should only perform a computation the first time it is called, and subsequent repeated calls are no-ops. If the user's intent is to re-compute the results (e.g., if parameters or other inputs have been changed), then setting this flag to true will force a new result to be computed.

Implements Isorropia::Operator.

Implemented in Isorropia::Epetra::Partitioner, Isorropia::Epetra::Colorer, Isorropia::Epetra::Partitioner2D, Isorropia::Epetra::Orderer, and Isorropia::Epetra::LevelScheduler.

bool Isorropia::Epetra::Operator::alreadyComputed ( ) const
inlinevirtual

Query whether compute_operation() has already been called.

Implements Isorropia::Operator.

int Isorropia::Epetra::Operator::numProperties ( ) const
inlinevirtual

Return the number of different values used for "properties".

For example, the number of colors or the number of parts used for the overall graph/matrix.

Returns
Global number of values for properties
Remarks
Infact, it returns the upper bound of the interval of taken values. For example, for the colors "1,2,4" , it will return "4"

Implements Isorropia::Operator.

int Isorropia::Epetra::Operator::numLocalProperties ( ) const
inlinevirtual

Return the number of different values used for "properties" for this process only.

Returns
Local number of values for properties

Implements Isorropia::Operator.

virtual const int& Isorropia::Epetra::Operator::operator[] ( int  myElem) const
virtual

Return the new partition ID for a given element that resided locally in the old operation.

Implements Isorropia::Operator.

virtual int Isorropia::Epetra::Operator::numElemsWithProperty ( int  property) const
virtual

Return the number of elements in a given partition.

Implements Isorropia::Operator.

virtual void Isorropia::Epetra::Operator::elemsWithProperty ( int  property,
int *  elementList,
int  len 
) const
virtual

Fill user-allocated list (of length len) with the global element ids to be located in the given partition.

Implements Isorropia::Operator.

virtual int Isorropia::Epetra::Operator::extractPropertiesCopy ( int  len,
int &  size,
int *  array 
) const
virtual

Copy a part of the property array.

Parameters
len[in] of the array given by the user.
size[out] Number of elements in the array.
array[out] Array of properties. Allocated by the user with a size of at least len elements.
Remarks
Memory space which is not useful in the array is not initialized or used in this method.
See Also
Isorropia::Operator::extractPropertiesView()

Implements Isorropia::Operator.

virtual int Isorropia::Epetra::Operator::extractPropertiesView ( int &  size,
const int *&  array 
) const
virtual

Give access of the property array that is owned by the current processor.

Parameters
size[out] Number of elements in the array.
array[out] Pointer to the the properties array inside the object.
Remarks
This pointer is only significant if the object still exists. Otherwise, you must use
See Also
Isorropia::Operator::extractPropertiesCopy().
Isorropia::Operator::extractPropertiesCopy()

Implements Isorropia::Operator.

void Isorropia::Epetra::Operator::paramsToUpper ( Teuchos::ParameterList &  ,
int &  changed,
bool  rmUnderscore = true 
)
private
void Isorropia::Epetra::Operator::stringToUpper ( std::string &  s,
int &  changed,
bool  rmUnderscore = false 
)
private
void Isorropia::Epetra::Operator::computeNumberOfProperties ( )
protected

Member Data Documentation

int Isorropia::Epetra::Operator::numberOfProperties_
private
int Isorropia::Epetra::Operator::localNumberOfProperties_
private
std::vector<int> Isorropia::Epetra::Operator::numberElemsByProperties_
private
Teuchos::RCP<const Epetra_BlockMap> Isorropia::Epetra::Operator::input_map_
protected
Teuchos::RCP<const Epetra_CrsGraph> Isorropia::Epetra::Operator::input_graph_
protected
Teuchos::RCP<const Epetra_RowMatrix> Isorropia::Epetra::Operator::input_matrix_
protected
Teuchos::RCP<const Epetra_MultiVector> Isorropia::Epetra::Operator::input_coords_
protected
Teuchos::RCP<Isorropia::Epetra::CostDescriber> Isorropia::Epetra::Operator::costs_
protected
Teuchos::RCP<const Epetra_MultiVector> Isorropia::Epetra::Operator::weights_
protected
Teuchos::ParameterList Isorropia::Epetra::Operator::paramlist_
protected
int Isorropia::Epetra::Operator::exportsSize_
protected
std::vector<int> Isorropia::Epetra::Operator::imports_
protected
std::vector<int> Isorropia::Epetra::Operator::properties_
protected
bool Isorropia::Epetra::Operator::operation_already_computed_
protected
int Isorropia::Epetra::Operator::global_num_vertex_weights_
protected
int Isorropia::Epetra::Operator::global_num_graph_edge_weights_
protected
int Isorropia::Epetra::Operator::global_num_hg_edge_weights_
protected
Teuchos::RCP<Library> Isorropia::Epetra::Operator::lib_
protected
int Isorropia::Epetra::Operator::base_
protected

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