196 TEUCHOS_TEST_FOR_EXCEPTION(
197 NumLocalParts_ < 1 || OverlappingLevel_ < 0,
199 "Ifpack2::OverlappingPartitioner::compute: "
200 "Invalid NumLocalParts_ or OverlappingLevel_.");
204 const char printMsg[] =
"OverlappingPartitioner: ";
206 if (verbose_ && (Graph_->getComm()->getRank() == 0)) {
207 cout << printMsg <<
"Number of local parts = "
208 << NumLocalParts_ << endl;
209 cout << printMsg <<
"Approx. Number of global parts = "
210 << NumLocalParts_ * Graph_->getComm ()->getSize () << endl;
211 cout << printMsg <<
"Amount of overlap = "
212 << OverlappingLevel_ << endl;
216 Partition_.resize (Graph_->getNodeNumRows ());
220 TEUCHOS_TEST_FOR_EXCEPTION(
221 ! Graph_->isFillComplete (), std::runtime_error,
222 "Ifpack2::OverlappingPartitioner::compute: "
223 "The input graph must be fill complete.");
225 TEUCHOS_TEST_FOR_EXCEPTION(
226 Graph_->getGlobalNumRows () != Graph_->getGlobalNumCols (),
228 "Ifpack2::OverlappingPartitioner::compute: "
229 "The input graph must be (globally) square.");
232 computePartitions ();
235 computeOverlappingPartitions ();
247 if (Partition_.size() == 0)
250 const local_ordinal_type invalid =
251 Teuchos::OrdinalTraits<local_ordinal_type>::invalid();
257 std::vector<size_t> sizes;
258 sizes.resize (NumLocalParts_);
261 for (
int i = 0; i < NumLocalParts_; ++i) {
265 for (
size_t i = 0; i < Graph_->getNodeNumRows (); ++i) {
266 TEUCHOS_TEST_FOR_EXCEPTION(
267 Partition_[i] >= NumLocalParts_, std::runtime_error,
268 "Ifpack2::OverlappingPartitioner::computeOverlappingPartitions: "
269 "Partition_[i] > NumLocalParts_.");
272 if (Partition_[i] != invalid) {
273 sizes[Partition_[i]]++;
278 Parts_.resize (NumLocalParts_);
279 for (
int i = 0; i < NumLocalParts_; ++i) {
280 Parts_[i].resize (sizes[i]);
284 for (
int i = 0; i < NumLocalParts_; ++i) {
288 for (
size_t i = 0; i < Graph_->getNodeNumRows (); ++i) {
289 const local_ordinal_type part = Partition_[i];
290 if (part != invalid) {
291 const size_t count = sizes[part];
292 Parts_[part][count] = i;
298 if (OverlappingLevel_ == 0) {
303 for (
int level = 1; level <= OverlappingLevel_; ++level) {
304 std::vector<std::vector<size_t> > tmp;
305 tmp.resize (NumLocalParts_);
311 int MaxNumEntries_tmp = Graph_->getNodeMaxNumRowEntries();
312 nonconst_local_inds_host_view_type Indices(
"Indices",MaxNumEntries_tmp);
313 nonconst_local_inds_host_view_type newIndices(
"newIndices",MaxNumEntries_tmp);
315 if (!maintainSparsity_) {
317 local_ordinal_type numLocalRows = Graph_->getNodeNumRows();
318 for (
int part = 0; part < NumLocalParts_ ; ++part) {
319 for (
size_t i = 0; i < Teuchos::as<size_t> (Parts_[part].size ()); ++i) {
320 const local_ordinal_type LRID = Parts_[part][i];
323 Graph_->getLocalRowCopy (LRID, Indices, numIndices);
325 for (
size_t j = 0; j < numIndices; ++j) {
327 const local_ordinal_type col = Indices[j];
328 if (col >= numLocalRows) {
333 std::vector<size_t>::iterator where =
334 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (col));
337 if (where == tmp[part].end()) {
338 tmp[part].push_back (col);
348 std::vector<size_t>::iterator where =
349 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (LRID));
353 if (where == tmp[part].end ()) {
354 tmp[part].push_back (LRID);
362 for (
int part = 0; part < NumLocalParts_ ; ++part) {
363 for (
size_t i = 0; i < Teuchos::as<size_t> (Parts_[part].size ()); ++i) {
364 const local_ordinal_type LRID = Parts_[part][i];
367 Graph_->getLocalRowCopy (LRID, Indices, numIndices);
372 Tpetra::sort(Indices,numIndices);
374 for (
size_t j = 0; j < numIndices; ++j) {
376 const local_ordinal_type col = Indices[j];
377 if (Teuchos::as<size_t> (col) >= Graph_->getNodeNumRows ()) {
382 std::vector<size_t>::iterator where =
383 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (col));
386 if (where == tmp[part].end()) {
389 size_t numNewIndices;
390 Graph_->getLocalRowCopy(col, newIndices, numNewIndices);
391 Tpetra::sort(newIndices,numNewIndices);
392 auto Indices_rcp = Kokkos::Compat::persistingView<nonconst_local_inds_host_view_type>(Indices, 0, numIndices);
393 auto newIndices_rcp = Kokkos::Compat::persistingView<nonconst_local_inds_host_view_type>(newIndices, 0, numNewIndices);
394 bool isSubset = std::includes(Indices_rcp.begin(),Indices_rcp.begin()+numIndices,
395 newIndices_rcp.begin(),newIndices_rcp.begin()+numNewIndices);
397 tmp[part].push_back (col);
403 std::vector<size_t>::iterator where =
404 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (LRID));
408 if (where == tmp[part].end ()) {
409 tmp[part].push_back (LRID);
421 for (
int i = 0; i < NumLocalParts_; ++i) {
422 Parts_[i].resize (tmp[i].size ());
423 for (
size_t j = 0; j < tmp[i].size (); ++j) {
424 Parts_[i][j] = tmp[i][j];