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

An implementation of the LevelScheduler interface that operates on and Epetra_CrsGraph, representing the non-zeros in a matrix. More...

#include <Isorropia_EpetraLevelScheduler.hpp>

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

Public Member Functions

 LevelScheduler (Teuchos::RCP< const Epetra_CrsGraph > input_graph, const Teuchos::ParameterList &paramlist=Teuchos::ParameterList("EmptyParameterList"), bool compute_now=true)
 Constructor. More...
 
 ~LevelScheduler ()
 Destructor. More...
 
void schedule (bool force_scheduling=false)
 Compute the scheduling if it has not already been computed, same effect as Isorropia::Epetra::LevelScheduler::compute. More...
 
void compute (bool force_compute=false)
 Compute the scheduling if it has not already been computed, same effect as Isorropia::Epetra::LevelScheduler::schedule. More...
 
virtual int numLevels () const
 Method which returns the number of levels. More...
 
virtual int numElemsWithLevel (int level) const
 Return the number of elements in a given level. More...
 
virtual void elemsWithLevel (int level, int *elementList, int len) const
 Fill user-allocated list (of length len) with the local ID for each element in the given level. More...
 
virtual void setParameters (const Teuchos::ParameterList &paramlist)=0
 Set parameters for the Operator instance. 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...
 
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...
 
virtual bool alreadyComputed () const =0
 Query whether the computation has already been called. More...
 
virtual int numProperties () const =0
 Return the number of different values used for "properties". More...
 
virtual int numLocalProperties () const =0
 Return the number of different values used for "properties" for this process only. More...
 
virtual const int & operator[] (int myElem) const =0
 Return the "property" for a given element that resided locally. More...
 
virtual int numElemsWithProperty (int property) const =0
 Return the number of LOCAL elements with the given property. More...
 
virtual void elemsWithProperty (int property, int *elementList, int len) const =0
 Fill user-allocated list (of length len) with the local element ids of the LOCAL elements with the given property. More...
 
virtual int extractPropertiesView (int &size, const int *&array) const =0
 Give access of the property array that is owned by the current processor. More...
 
virtual int extractPropertiesCopy (int len, int &size, int *array) const =0
 Copy a part of the property array. 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_
 

Detailed Description

An implementation of the LevelScheduler interface that operates on and Epetra_CrsGraph, representing the non-zeros in a matrix.

The elements to be partitioned into levels are matrix rows. Assumption is that matrix is lower triangular or upper triangular.

Constructor & Destructor Documentation

Isorropia::Epetra::LevelScheduler::LevelScheduler ( Teuchos::RCP< const Epetra_CrsGraph >  input_graph,
const Teuchos::ParameterList &  paramlist = Teuchos::ParameterList("EmptyParameterList"),
bool  compute_now = true 
)

Constructor.

Parameters
input_graph(in) the graph representing the non-zeros of the matrix
paramlist(in) list of parameters
compute_now(in) if true, the scheduling is computed in the constructor, otherwise call Isorropia::Epetra::LevelScheduler::schedule when you want to compute the scheduling, defaults to false
Isorropia::Epetra::LevelScheduler::~LevelScheduler ( )
virtual

Destructor.

Reimplemented from Isorropia::LevelScheduler.

Member Function Documentation

void Isorropia::Epetra::LevelScheduler::schedule ( bool  force_scheduling = false)
virtual

Compute the scheduling if it has not already been computed, same effect as Isorropia::Epetra::LevelScheduler::compute.

Parameters
force_scheduling(in) if true recompute the scheduling even if it has already been computed, defaults to false

Implements Isorropia::LevelScheduler.

void Isorropia::Epetra::LevelScheduler::compute ( bool  force_compute = false)
inlinevirtual

Compute the scheduling if it has not already been computed, same effect as Isorropia::Epetra::LevelScheduler::schedule.

Parameters
force_compute(in) if true recompute the scheduling even if it has already been computed, defaults to false

Implements Isorropia::Epetra::Operator.

virtual int Isorropia::LevelScheduler::numLevels ( ) const
inlinevirtualinherited

Method which returns the number of levels.

Returns
The number of levels on the local process. Levels begin at 1 and increase by 1.
See Also
Isorropia::Operator::numProperties()
virtual int Isorropia::LevelScheduler::numElemsWithLevel ( int  level) const
inlinevirtualinherited

Return the number of elements in a given level.

Parameters
levelThe wanted level.
Returns
The number of elements in the level.
See Also
Isorropia::Operator::numElemsWithProperty()
virtual void Isorropia::LevelScheduler::elemsWithLevel ( int  level,
int *  elementList,
int  len 
) const
inlinevirtualinherited

Fill user-allocated list (of length len) with the local ID for each element in the given level.

Parameters
levelthe wanted level
elementListan array to receive local elements of the given level
lenthe number of elements wanted
See Also
Isorropia::Operator::elemsWithProperty()
virtual void Isorropia::Operator::setParameters ( const Teuchos::ParameterList &  paramlist)
pure virtualinherited

Set parameters for the Operator 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.

Parameters
paramlist[in] List of parameters that the user wants to use.

Implemented in Isorropia::Epetra::Operator.

virtual bool Isorropia::Operator::alreadyComputed ( ) const
pure virtualinherited

Query whether the computation has already been called.

Returns
True if the computation has already been done, False otherwise.

Implemented in Isorropia::Epetra::Operator.

virtual int Isorropia::Operator::numProperties ( ) const
pure virtualinherited

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"

Implemented in Isorropia::Epetra::Operator.

virtual int Isorropia::Operator::numLocalProperties ( ) const
pure virtualinherited

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

Returns
Local number of values for properties

Implemented in Isorropia::Epetra::Operator.

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

Return the "property" for a given element that resided locally.

Parameters
myElem[in] the local ID of the element we want to know the property.
Returns
property associated to the local element.

Implemented in Isorropia::Epetra::Operator.

virtual int Isorropia::Operator::numElemsWithProperty ( int  property) const
pure virtualinherited

Return the number of LOCAL elements with the given property.

Parameters
property[in] Value of the property to consider.
Returns
Number of local elems which have this property.

Implemented in Isorropia::Epetra::Operator.

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

Fill user-allocated list (of length len) with the local element ids of the LOCAL elements with the given property.

Parameters
property[in] Value of the property to consider.
elementList[out] User allocated array (of size at least len) of local ID that have the asked property.
len[in] Maximum lenght for the array. If len is greater than the result of numElemsWithProperty() for property, only the first and relevant elements are filled.
Remarks
Memory space which is not useful in the array is not initialized or used in this method.

Implemented in Isorropia::Epetra::Operator.

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

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()

Implemented in Isorropia::Epetra::Operator.

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

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()

Implemented in Isorropia::Epetra::Operator.

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

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 ( )
inlineinherited

Get the cost object.

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

Query whether compute_operation() has already been called.

Implements Isorropia::Operator.

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

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
inlinevirtualinherited

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
virtualinherited

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
virtualinherited

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
virtualinherited

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
virtualinherited

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
virtualinherited

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::computeNumberOfProperties ( )
protectedinherited

Member Data Documentation

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

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