Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_GreedyPartitioner.cpp
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #include "Ifpack_ConfigDefs.h"
44 #include "Ifpack_Partitioner.h"
47 #include "Ifpack_Graph.h"
48 
49 #include "Epetra_Comm.h"
50 #include "Epetra_BlockMap.h"
51 #include "Epetra_Map.h"
52 #include "Epetra_CrsGraph.h"
53 #include "Teuchos_ParameterList.hpp"
54 
55 //==============================================================================
57 {
58  std::vector<int> ElementsPerPart(NumLocalParts());
59  std::vector<int> Count(NumLocalParts());
60  for (int i = 0 ; i < NumLocalParts() ; ++i)
61  Count[i] = 0;
62 
63  // define how many nodes have to be put on each part
64  int div = NumMyRows() / NumLocalParts();
65  int mod = NumMyRows() % NumLocalParts();
66 
67  for (int i = 0 ; i < NumLocalParts() ; ++i) {
68  Count[i] = 0;
69  ElementsPerPart[i] = div;
70  if (i < mod) ElementsPerPart[i]++;
71  }
72 
73  for( int i=0 ; i<NumMyRows() ; ++i ) {
74  Partition_[i] = -1;
75  }
76 
77  int NumEntries;
78  std::vector<int> Indices(MaxNumEntries());
79 
80  // load root node for partition 0
81  int CurrentPartition = 0;
82  int TotalCount = 0;
83 
84  // filter singletons and empty rows, put all of them in partition 0
85  for (int i = 0 ; i < NumMyRows() ; ++i) {
86  NumEntries = 0;
87  int ierr = Graph_->ExtractMyRowCopy(i, MaxNumEntries(),
88  NumEntries, &Indices[0]);
89  IFPACK_CHK_ERR(ierr);
90  if (NumEntries <= 1) {
91  Partition_[i] = 0;
92  TotalCount++;
93  }
94  }
95 
96  if (TotalCount)
97  CurrentPartition = 1;
98 
99  std::vector<int> ThisLevel(1);
100  ThisLevel[0] = RootNode_;
101 
102  // be sure that RootNode is not a singleton or empty row
103  if (Partition_[RootNode_] != -1) {
104  // look for another RN
105  for (int i = 0 ; i < NumMyRows() ; ++i)
106  if (Partition_[i] == -1) {
107  ThisLevel[0] = i;
108  break;
109  }
110  }
111  else {
112  Partition_[RootNode_] = CurrentPartition;
113  }
114 
115  // now aggregate the non-empty and non-singleton rows
116  while (ThisLevel.size()) {
117 
118  std::vector<int> NextLevel;
119 
120  for (unsigned int i = 0 ; i < ThisLevel.size() ; ++i) {
121 
122  int CurrentNode = ThisLevel[i];
123  int ierr = Graph_->ExtractMyRowCopy(CurrentNode, MaxNumEntries(),
124  NumEntries, &Indices[0]);
125  IFPACK_CHK_ERR(ierr);
126 
127  if (NumEntries <= 1)
128  continue;
129 
130  for (int j = 0 ; j < NumEntries ; ++j) {
131 
132  int NextNode = Indices[j];
133  if (NextNode >= NumMyRows()) continue;
134 
135  if (Partition_[NextNode] == -1) {
136  // this is a free node
137  NumLocalParts_ = CurrentPartition + 1;
138  Partition_[NextNode] = CurrentPartition;
139  ++Count[CurrentPartition];
140  ++TotalCount;
141  NextLevel.push_back(NextNode);
142  }
143  }
144  } // for (i)
145 
146  // check whether change partition or not
147  if (Count[CurrentPartition] >= ElementsPerPart[CurrentPartition])
148  ++CurrentPartition;
149 
150  // swap next and this
151  ThisLevel.resize(0);
152  for (unsigned int i = 0 ; i < NextLevel.size() ; ++i)
153  ThisLevel.push_back(NextLevel[i]);
154 
155  if (ThisLevel.size() == 0 && (TotalCount != NumMyRows())) {
156  // need to look for new RootNode, do this in a simple way
157  for (int i = 0 ; i < NumMyRows() ; i++) {
158  if (Partition_[i] == -1)
159  ThisLevel.push_back(i);
160  break;
161  }
162  }
163 
164  } // while (ok)
165 
166  return(0);
167 }
168 
virtual int ExtractMyRowCopy(int MyRow, int LenOfIndices, int &NumIndices, int *Indices) const =0
Extracts a copy of input local row.
int MaxNumEntries() const
Returns the max number of local entries in a row.
int ComputePartitions()
Computes the partitions. Returns 0 if successful.
int NumMyRows() const
Returns the number of local rows.
std::vector< int > Partition_
Partition_[i] contains the ID of non-overlapping part it belongs to.
const Ifpack_Graph * Graph_
Reference to the graph to be partitioned.
int NumLocalParts_
Number of local subgraphs.
int NumLocalParts() const
Returns the number of computed local partitions.
#define IFPACK_CHK_ERR(ifpack_err)