242 using Teuchos::Array;
243 using Teuchos::ArrayView;
246 using Teuchos::REDUCE_SUM;
247 using Teuchos::reduceAll;
249 size_t NumIn, NumL, NumU;
252 constructOverlapGraph();
255 const int MaxNumIndices = OverlapGraph_->getLocalMaxNumRowEntries ();
259 const int NumMyRows = OverlapGraph_->getRowMap ()->getLocalNumElements ();
261 using device_type =
typename node_type::device_type;
262 using execution_space =
typename device_type::execution_space;
263 using dual_view_type = Kokkos::DualView<size_t*,device_type>;
264 dual_view_type numEntPerRow_dv(
"numEntPerRow",NumMyRows);
265 Tpetra::Details::WrappedDualView<dual_view_type> numEntPerRow(numEntPerRow_dv);
267 const auto overalloc = Overalloc_;
268 const auto levelfill = LevelFill_;
271 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
272 auto localOverlapGraph = OverlapGraph_->getLocalGraphDevice();
273 Kokkos::parallel_for(
"CountOverlapGraphRowEntries",
274 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
275 KOKKOS_LAMBDA(
const int i)
278 int RowMaxNumIndices = localOverlapGraph.rowConst(i).length;
279 numEntPerRow_d(i) = (levelfill == 0) ? RowMaxNumIndices
280 : Kokkos::ceil(
static_cast<double>(RowMaxNumIndices)
281 * Kokkos::pow(overalloc, levelfill));
289 Teuchos::ArrayView<const size_t> a_numEntPerRow(numEntPerRow.getHostView(Tpetra::Access::ReadOnly).data(),NumMyRows);
291 OverlapGraph_->getRowMap (),
294 OverlapGraph_->getRowMap (),
297 Array<local_ordinal_type> L (MaxNumIndices);
298 Array<local_ordinal_type> U (MaxNumIndices);
304 for (
int i = 0; i< NumMyRows; ++i) {
305 local_inds_host_view_type my_indices;
306 OverlapGraph_->getLocalRowView (i, my_indices);
313 NumIn = my_indices.size();
315 for (
size_t j = 0; j < NumIn; ++j) {
316 const local_ordinal_type k = my_indices[j];
340 L_Graph_->insertLocalIndices (i, NumL, L.data());
343 U_Graph_->insertLocalIndices (i, NumU, U.data());
347 if (LevelFill_ > 0) {
349 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
350 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
351 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
352 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
353 RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
354 params->set (
"Optimize Storage",
false);
355 L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
356 U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
357 L_Graph_->resumeFill ();
358 U_Graph_->resumeFill ();
364 int MaxRC = NumMyRows;
365 std::vector<std::vector<int> > Levels(MaxRC);
366 std::vector<int> LinkList(MaxRC);
367 std::vector<int> CurrentLevel(MaxRC);
368 Array<local_ordinal_type> CurrentRow (MaxRC + 1);
369 std::vector<int> LevelsRowU(MaxRC);
372 for (
int i = 0; i < NumMyRows; ++i) {
377 size_t LenL = L_Graph_->getNumEntriesInLocalRow(i);
378 size_t LenU = U_Graph_->getNumEntriesInLocalRow(i);
379 size_t Len = LenL + LenU + 1;
380 CurrentRow.resize(Len);
381 nonconst_local_inds_host_view_type CurrentRow_view(CurrentRow.data(),CurrentRow.size());
382 L_Graph_->getLocalRowCopy(i, CurrentRow_view, LenL);
383 CurrentRow[LenL] = i;
385 ArrayView<local_ordinal_type> URowView = CurrentRow.view (LenL+1,LenU);
386 nonconst_local_inds_host_view_type URowView_v(URowView.data(),URowView.size());
389 U_Graph_->getLocalRowCopy (i, URowView_v, LenU);
394 for (
size_t j=0; j<Len-1; j++) {
395 LinkList[CurrentRow[j]] = CurrentRow[j+1];
396 CurrentLevel[CurrentRow[j]] = 0;
399 LinkList[CurrentRow[Len-1]] = NumMyRows;
400 CurrentLevel[CurrentRow[Len-1]] = 0;
404 First = CurrentRow[0];
407 int PrevInList = Next;
408 int NextInList = LinkList[Next];
411 local_inds_host_view_type IndicesU;
412 U_Graph_->getLocalRowView (RowU, IndicesU);
414 int LengthRowU = IndicesU.size ();
420 for (ii = 0; ii < LengthRowU; ) {
421 int CurInList = IndicesU[ii];
422 if (CurInList < NextInList) {
424 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
425 if (NewLevel <= LevelFill_) {
426 LinkList[PrevInList] = CurInList;
427 LinkList[CurInList] = NextInList;
428 PrevInList = CurInList;
429 CurrentLevel[CurInList] = NewLevel;
433 else if (CurInList == NextInList) {
434 PrevInList = NextInList;
435 NextInList = LinkList[PrevInList];
436 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
437 CurrentLevel[CurInList] = std::min (CurrentLevel[CurInList],
442 PrevInList = NextInList;
443 NextInList = LinkList[PrevInList];
446 Next = LinkList[Next];
450 CurrentRow.resize(0);
456 CurrentRow.push_back(Next);
457 Next = LinkList[Next];
463 L_Graph_->removeLocalIndices (i);
464 if (CurrentRow.size() > 0) {
465 L_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
470 TEUCHOS_TEST_FOR_EXCEPTION(
471 Next != i, std::runtime_error,
472 "Ifpack2::IlukGraph::initialize: FATAL: U has zero diagonal")
474 LevelsRowU[0] = CurrentLevel[Next];
475 Next = LinkList[Next];
478 CurrentRow.resize(0);
481 while (Next < NumMyRows) {
482 LevelsRowU[LenU+1] = CurrentLevel[Next];
483 CurrentRow.push_back (Next);
485 Next = LinkList[Next];
492 U_Graph_->removeLocalIndices (i);
494 U_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
498 Levels[i] = std::vector<int> (LenU+1);
499 for (
size_t jj=0; jj<LenU+1; jj++) {
500 Levels[i][jj] = LevelsRowU[jj];
504 catch (std::runtime_error &e) {
506 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
507 Kokkos::parallel_for(
"CountOverlapGraphRowEntries",
508 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
509 KOKKOS_LAMBDA(
const int i)
511 const auto numRowEnt = numEntPerRow_d(i);
512 numEntPerRow_d(i) = ceil(
static_cast<double>((numRowEnt != 0 ? numRowEnt : 1)) * overalloc);
515 const int localInsertError = insertError ? 1 : 0;
516 int globalInsertError = 0;
517 reduceAll (* (OverlapGraph_->getRowMap ()->getComm ()), REDUCE_SUM, 1,
518 &localInsertError, &globalInsertError);
519 insertError = globalInsertError > 0;
521 }
while (insertError);
524 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
525 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
526 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
527 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
528 L_Graph_->fillComplete (L_DomainMap, L_RangeMap);
529 U_Graph_->fillComplete (U_DomainMap, U_RangeMap);
531 reduceAll<int, size_t> (* (L_DomainMap->getComm ()), REDUCE_SUM, 1,
532 &NumMyDiagonals_, &NumGlobalDiagonals_);