129 NumLocalParts_ = List.get(
"partitioner: local parts", NumLocalParts_);
130 OverlappingLevel_ = List.get(
"partitioner: overlap", OverlappingLevel_);
131 verbose_ = List.get(
"partitioner: print level", verbose_);
132 maintainSparsity_ = List.get(
"partitioner: maintain sparsity",
false);
133 typedef Teuchos::RCP< Tpetra::Map<typename GraphType::local_ordinal_type, typename GraphType::global_ordinal_type, typename GraphType::node_type>
const > map_type;
134 typedef Teuchos::RCP<Tpetra::Import<typename GraphType::local_ordinal_type, typename GraphType::global_ordinal_type, typename GraphType::node_type>
const > import_type;
139 import_type theImport = Graph_->getImporter();
140 if (theImport != Teuchos::null) List.set< import_type >(
"theImport",theImport);
141 List.set< map_type >(
"OverlapRowMap",Graph_->getRowMap());
143 if (NumLocalParts_ < 0) {
144 NumLocalParts_ = Graph_->getLocalNumRows() / (-NumLocalParts_);
146 if (NumLocalParts_ == 0) {
151 TEUCHOS_TEST_FOR_EXCEPTION(
152 NumLocalParts_ < 0 ||
153 Teuchos::as<size_t> (NumLocalParts_) > Graph_->getLocalNumRows(),
155 "Ifpack2::OverlappingPartitioner::setParameters: "
156 "Invalid NumLocalParts_ = " << NumLocalParts_ <<
".");
157 TEUCHOS_TEST_FOR_EXCEPTION(
158 OverlappingLevel_ < 0, std::runtime_error,
159 "Ifpack2::OverlappingPartitioner::setParameters: "
160 "Invalid OverlappingLevel_ = " << OverlappingLevel_ <<
".");
162 setPartitionParameters(List);
172 TEUCHOS_TEST_FOR_EXCEPTION(
173 NumLocalParts_ < 1 || OverlappingLevel_ < 0,
175 "Ifpack2::OverlappingPartitioner::compute: "
176 "Invalid NumLocalParts_ or OverlappingLevel_.");
180 const char printMsg[] =
"OverlappingPartitioner: ";
182 if (verbose_ && (Graph_->getComm()->getRank() == 0)) {
183 cout << printMsg <<
"Number of local parts = "
184 << NumLocalParts_ << endl;
185 cout << printMsg <<
"Approx. Number of global parts = "
186 << NumLocalParts_ * Graph_->getComm ()->getSize () << endl;
187 cout << printMsg <<
"Amount of overlap = "
188 << OverlappingLevel_ << endl;
192 Partition_.resize (Graph_->getLocalNumRows ());
196 TEUCHOS_TEST_FOR_EXCEPTION(
197 ! Graph_->isFillComplete (), std::runtime_error,
198 "Ifpack2::OverlappingPartitioner::compute: "
199 "The input graph must be fill complete.");
201 TEUCHOS_TEST_FOR_EXCEPTION(
202 Graph_->getGlobalNumRows () != Graph_->getGlobalNumCols (),
204 "Ifpack2::OverlappingPartitioner::compute: "
205 "The input graph must be (globally) square.");
208 computePartitions ();
211 computeOverlappingPartitions ();
223 if (Partition_.size() == 0)
226 const local_ordinal_type invalid =
227 Teuchos::OrdinalTraits<local_ordinal_type>::invalid();
233 std::vector<size_t> sizes;
234 sizes.resize (NumLocalParts_);
237 for (
int i = 0; i < NumLocalParts_; ++i) {
241 for (
size_t i = 0; i < Graph_->getLocalNumRows (); ++i) {
242 TEUCHOS_TEST_FOR_EXCEPTION(
243 Partition_[i] >= NumLocalParts_, std::runtime_error,
244 "Ifpack2::OverlappingPartitioner::computeOverlappingPartitions: "
245 "Partition_[i] > NumLocalParts_.");
248 if (Partition_[i] != invalid) {
249 sizes[Partition_[i]]++;
254 Parts_.resize (NumLocalParts_);
255 for (
int i = 0; i < NumLocalParts_; ++i) {
256 Parts_[i].resize (sizes[i]);
260 for (
int i = 0; i < NumLocalParts_; ++i) {
264 for (
size_t i = 0; i < Graph_->getLocalNumRows (); ++i) {
265 const local_ordinal_type part = Partition_[i];
266 if (part != invalid) {
267 const size_t count = sizes[part];
268 Parts_[part][count] = i;
274 if (OverlappingLevel_ == 0) {
279 for (
int level = 1; level <= OverlappingLevel_; ++level) {
280 std::vector<std::vector<size_t> > tmp;
281 tmp.resize (NumLocalParts_);
287 int MaxNumEntries_tmp = Graph_->getLocalMaxNumRowEntries();
288 nonconst_local_inds_host_view_type Indices(
"Indices",MaxNumEntries_tmp);
289 nonconst_local_inds_host_view_type newIndices(
"newIndices",MaxNumEntries_tmp);
291 if (!maintainSparsity_) {
293 local_ordinal_type numLocalRows = Graph_->getLocalNumRows();
294 for (
int part = 0; part < NumLocalParts_ ; ++part) {
295 for (
size_t i = 0; i < Teuchos::as<size_t> (Parts_[part].size ()); ++i) {
296 const local_ordinal_type LRID = Parts_[part][i];
299 Graph_->getLocalRowCopy (LRID, Indices, numIndices);
301 for (
size_t j = 0; j < numIndices; ++j) {
303 const local_ordinal_type col = Indices[j];
304 if (col >= numLocalRows) {
309 std::vector<size_t>::iterator where =
310 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (col));
313 if (where == tmp[part].end()) {
314 tmp[part].push_back (col);
324 std::vector<size_t>::iterator where =
325 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (LRID));
329 if (where == tmp[part].end ()) {
330 tmp[part].push_back (LRID);
338 for (
int part = 0; part < NumLocalParts_ ; ++part) {
339 for (
size_t i = 0; i < Teuchos::as<size_t> (Parts_[part].size ()); ++i) {
340 const local_ordinal_type LRID = Parts_[part][i];
343 Graph_->getLocalRowCopy (LRID, Indices, numIndices);
348 Tpetra::sort(Indices,numIndices);
350 for (
size_t j = 0; j < numIndices; ++j) {
352 const local_ordinal_type col = Indices[j];
353 if (Teuchos::as<size_t> (col) >= Graph_->getLocalNumRows ()) {
358 std::vector<size_t>::iterator where =
359 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (col));
362 if (where == tmp[part].end()) {
365 size_t numNewIndices;
366 Graph_->getLocalRowCopy(col, newIndices, numNewIndices);
367 Tpetra::sort(newIndices,numNewIndices);
368 auto Indices_rcp = Kokkos::Compat::persistingView<nonconst_local_inds_host_view_type>(Indices, 0, numIndices);
369 auto newIndices_rcp = Kokkos::Compat::persistingView<nonconst_local_inds_host_view_type>(newIndices, 0, numNewIndices);
370 bool isSubset = std::includes(Indices_rcp.begin(),Indices_rcp.begin()+numIndices,
371 newIndices_rcp.begin(),newIndices_rcp.begin()+numNewIndices);
373 tmp[part].push_back (col);
379 std::vector<size_t>::iterator where =
380 std::find (tmp[part].begin (), tmp[part].end (), Teuchos::as<size_t> (LRID));
384 if (where == tmp[part].end ()) {
385 tmp[part].push_back (LRID);
397 for (
int i = 0; i < NumLocalParts_; ++i) {
398 Parts_[i].resize (tmp[i].size ());
399 for (
size_t j = 0; j < tmp[i].size (); ++j) {
400 Parts_[i][j] = tmp[i][j];