00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "Ifpack_ConfigDefs.h"
00044 #include "Ifpack_Partitioner.h"
00045 #include "Ifpack_OverlappingPartitioner.h"
00046 #include "Ifpack_Graph.h"
00047
00048 #include "Epetra_Comm.h"
00049 #include "Epetra_BlockMap.h"
00050 #include "Epetra_Map.h"
00051 #include "Teuchos_ParameterList.hpp"
00052
00053 static const string PrintMsg_ = "(Ifpack_OvPartitioner) ";
00054
00055
00056 Ifpack_OverlappingPartitioner::
00057 Ifpack_OverlappingPartitioner(const Ifpack_Graph* Graph) :
00058 NumLocalParts_(1),
00059 Graph_(Graph),
00060 OverlappingLevel_(0),
00061 IsComputed_(false),
00062 verbose_(false)
00063 {
00064 }
00065
00066
00067 Ifpack_OverlappingPartitioner::~Ifpack_OverlappingPartitioner()
00068 {
00069 }
00070
00071
00072 int Ifpack_OverlappingPartitioner::SetParameters(Teuchos::ParameterList& List)
00073 {
00074
00075 NumLocalParts_ = List.get("partitioner: local parts", NumLocalParts_);
00076 OverlappingLevel_ = List.get("partitioner: overlap", OverlappingLevel_);
00077 verbose_ = List.get("partitioner: print level", verbose_);
00078
00079 if (NumLocalParts_ < 0)
00080 NumLocalParts_ = Graph_->NumMyRows() / (-NumLocalParts_);
00081 if (NumLocalParts_ == 0)
00082 NumLocalParts_ = 1;
00083 if (NumLocalParts_ < 0)
00084 IFPACK_CHK_ERR(-1);
00085 if (NumLocalParts_ > Graph_->NumMyRows())
00086 IFPACK_CHK_ERR(-1);
00087
00088 if (OverlappingLevel_ < 0)
00089 IFPACK_CHK_ERR(-1);
00090
00091 SetPartitionParameters(List);
00092
00093 return(0);
00094 }
00095
00096
00097 int Ifpack_OverlappingPartitioner::Compute()
00098 {
00099
00100 if (NumLocalParts_ < 1)
00101 IFPACK_CHK_ERR(-1);
00102
00103 if (OverlappingLevel_ < 0)
00104 IFPACK_CHK_ERR(-1);
00105
00106
00107
00108 if (verbose_ && (Comm().MyPID() == 0)) {
00109 cout << PrintMsg_ << "Number of local parts = " << NumLocalParts_ << endl;
00110 cout << PrintMsg_ << "Number of global parts = "
00111 << NumLocalParts_ * Comm().NumProc() << endl;
00112 cout << PrintMsg_ << "Amount of overlap = " << OverlappingLevel_ << endl;
00113 }
00114
00115
00116
00117 Partition_.resize(NumMyRows());
00118 Parts_.resize(NumLocalParts());
00119
00120
00121
00122 if (Graph_->Filled() == false)
00123 IFPACK_CHK_ERR(-4);
00124
00125 if (Graph_->NumGlobalRows64() != Graph_->NumGlobalCols64())
00126 IFPACK_CHK_ERR(-3);
00127
00128 if (NumLocalParts_ < 1)
00129 IFPACK_CHK_ERR(-2);
00130
00131
00132
00133 IFPACK_CHK_ERR(ComputePartitions());
00134
00135
00136
00137 IFPACK_CHK_ERR(ComputeOverlappingPartitions());
00138
00139
00140
00141 IsComputed_ = true;
00142
00143 return(0);
00144 }
00145
00146
00147 int Ifpack_OverlappingPartitioner::ComputeOverlappingPartitions()
00148 {
00149
00150
00151
00152
00153 std::vector<int> sizes;
00154 sizes.resize(NumLocalParts_);
00155
00156
00157 for (int i = 0 ; i < NumLocalParts_ ; ++i)
00158 sizes[i] = 0;
00159
00160 for (int i = 0 ; i < NumMyRows() ; ++i) {
00161 if (Partition_[i] >= NumLocalParts_) {
00162 cerr << "ERROR: Partition[" << i << "] = "<< Partition_[i]
00163 << ", NumLocalParts = " << NumLocalParts_ << endl;
00164 cerr << "(file = " << __FILE__ << ", line = "
00165 << __LINE__ << ")" << endl;
00166 IFPACK_CHK_ERR(-10);
00167 }
00168
00169
00170 if (Partition_[i] == -1)
00171 IFPACK_CHK_ERR(-1);
00172 sizes[Partition_[i]]++;
00173 }
00174
00175
00176 for (int i = 0 ; i < NumLocalParts_ ; ++i)
00177 Parts_[i].resize(sizes[i]);
00178
00179
00180
00181 for (int i = 0 ; i < NumLocalParts_ ; ++i)
00182 sizes[i] = 0;
00183
00184 for (int i = 0 ; i < NumMyRows() ; ++i) {
00185 int part = Partition_[i];
00186 int count = sizes[part];
00187 Parts_[part][count] = i;
00188 sizes[part]++;
00189 }
00190
00191 if (OverlappingLevel_ == 0)
00192 return(0);
00193
00194
00195 for (int level = 1 ; level <= OverlappingLevel_ ; ++level) {
00196
00197 std::vector<std::vector<int> > tmp;
00198 tmp.resize(NumLocalParts_);
00199
00200
00201
00202
00203
00204 int MaxNumEntries_tmp = Graph_->MaxMyNumEntries();
00205 std::vector<int> Indices;
00206 Indices.resize(MaxNumEntries_tmp);
00207
00208 for (int part = 0 ; part < NumLocalParts_ ; ++part) {
00209
00210 for (int i = 0; i < (int)Parts_[part].size() ; ++i) {
00211
00212 int LRID = Parts_[part][i];
00213 int NumIndices;
00214 int ierr = Graph_->ExtractMyRowCopy(LRID, MaxNumEntries_tmp,
00215 NumIndices, &Indices[0]);
00216 IFPACK_CHK_ERR(ierr);
00217
00218 for (int j = 0 ; j < NumIndices ; ++j) {
00219
00220
00221 int col = Indices[j];
00222 if (col >= NumMyRows())
00223 continue;
00224
00225
00226 std::vector<int>::iterator
00227 where = find(tmp[part].begin(), tmp[part].end(), col);
00228
00229 if (where == tmp[part].end()) {
00230 tmp[part].push_back(col);
00231 }
00232 }
00233 }
00234 }
00235
00236
00237 for (int i = 0 ; i < NumLocalParts_ ; ++i) {
00238 Parts_[i].resize(tmp[i].size());
00239 for (int j = 0 ; j < (int)tmp[i].size() ; ++j)
00240 Parts_[i][j] = tmp[i][j];
00241 }
00242 }
00243
00244 return(0);
00245
00246 }
00247
00248
00249 int Ifpack_OverlappingPartitioner::NumMyRows() const
00250 {
00251 return(Graph_->NumMyRows());
00252 }
00253
00254
00255 int Ifpack_OverlappingPartitioner::NumMyNonzeros() const
00256 {
00257 return(Graph_->NumMyNonzeros());
00258 }
00259
00260
00261 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00262 int Ifpack_OverlappingPartitioner::NumGlobalRows() const
00263 {
00264 return(Graph_->NumGlobalRows());
00265 }
00266 #endif
00267
00268 long long Ifpack_OverlappingPartitioner::NumGlobalRows64() const
00269 {
00270 return(Graph_->NumGlobalRows64());
00271 }
00272
00273 int Ifpack_OverlappingPartitioner::MaxNumEntries() const
00274 {
00275 return(Graph_->MaxMyNumEntries());
00276 }
00277
00278
00279 const Epetra_Comm& Ifpack_OverlappingPartitioner::Comm() const
00280 {
00281 return(Graph_->Comm());
00282 }
00283
00284
00285 ostream& Ifpack_OverlappingPartitioner::Print(ostream & os) const
00286 {
00287
00288 if (Comm().MyPID())
00289 return(os);
00290
00291 os << "================================================================================" << endl;
00292 os << "Ifpack_OverlappingPartitioner" << endl;
00293 os << "Number of local rows = " << Graph_->NumMyRows() << endl;
00294 os << "Number of global rows = " << Graph_->NumGlobalRows64() << endl;
00295 os << "Number of local parts = " << NumLocalParts_ << endl;
00296 os << "Overlapping level = " << OverlappingLevel_ << endl;
00297 os << "Is computed = " << IsComputed_ << endl;
00298 os << "================================================================================" << endl;
00299
00300 return(os);
00301 }
00302