43 #ifndef IFPACK_OVERLAPPINGPARTITIONER_H 44 #define IFPACK_OVERLAPPINGPARTITIONER_H 48 #include "Teuchos_ParameterList.hpp" 52 class Epetra_BlockMap;
113 if ((MyRow < 0) || (MyRow >
NumMyRows()))
125 if ((j < 0) || (j > (
int)
Parts_[i].size()))
134 return(
Parts_[Part].size());
140 List[i] =
Parts_[Part][i];
181 virtual std::ostream&
Print(std::ostream& os)
const;
189 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 197 const Epetra_Comm&
Comm()
const;
216 #endif // IFPACK_OVERLAPPINGPARTITIONER_H const int * NonOverlappingPartition() const
Returns a pointer to the integer vector containing the non-overlapping partition ID of each local row...
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the partitioner.
virtual ~Ifpack_OverlappingPartitioner()
Destructor.
virtual int Compute()
Computes the partitions. Returns 0 if successful.
bool IsComputed()
Returns true if partitions have been computed successfully.
int NumGlobalRows() const
Returns the number of global rows.
virtual int SetPartitionParameters(Teuchos::ParameterList &List)=0
Sets all the parameters for the partitioner.
Ifpack_OverlappingPartitioner(const Ifpack_Graph *Graph)
Constructor.
int RowsInPart(const int Part, int *List) const
Copies into List the rows in the (overlapping) partition Part.
std::vector< int > Partition_
Partition_[i] contains the ID of non-overlapping part it belongs to.
int OverlappingLevel_
Overlapping level.
Ifpack_Partitioner: A class to decompose local Ifpack_Graph's.
int operator()(int MyRow) const
Returns the local non-overlapping partition ID of the specified row.
virtual int ComputeOverlappingPartitions()
Computes the partitions. Returns 0 if successful.
bool IsComputed_
If true, the graph has been successfully partitioned.
long long NumGlobalRows64() const
bool verbose_
If true, information are reported on cout.
int NumMyRows() const
Returns the number of local rows.
int NumLocalParts() const
Returns the number of computed local partitions.
adjacency_list< vecS, vecS, undirectedS, no_property, property< edge_weight_t, double > > Graph
int OverlappingLevel() const
Returns the overlapping level.
const Ifpack_Graph * Graph_
Reference to the graph to be partitioned.
int NumRowsInPart(const int Part) const
Returns the number of rows contained in specified partition.
int NumLocalParts_
Number of local subgraphs.
virtual int ComputePartitions()=0
Computes the partitions. Returns 0 if successful.
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator<<.
Ifpack_Graph: a pure virtual class that defines graphs for IFPACK.
std::vector< std::vector< int > > Parts_
Parts_[i][j] is the ID of the j-th row contained in the (overlapping)
const Epetra_Comm & Comm() const
Returns the communicator object of Graph.
#define IFPACK_CHK_ERR(ifpack_err)
int MaxNumEntries() const
Returns the max number of local entries in a row.
int NumMyNonzeros() const
Returns the number of local nonzero elements.