Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_IlukGraph.hpp
Go to the documentation of this file.
1/*@HEADER
2// ***********************************************************************
3//
4// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5// Copyright (2009) 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
48
49#ifndef IFPACK2_ILUK_GRAPH_HPP
50#define IFPACK2_ILUK_GRAPH_HPP
51
52#include <algorithm>
53#include <vector>
54
55#include "KokkosSparse_spiluk.hpp"
56
57#include <Ifpack2_ConfigDefs.hpp>
58#include <Teuchos_ParameterList.hpp>
59#include <Teuchos_CommHelpers.hpp>
60#include <Tpetra_CrsGraph.hpp>
61#include <Tpetra_Details_WrappedDualView.hpp>
62#include <Tpetra_Import.hpp>
63#include <Ifpack2_CreateOverlapGraph.hpp>
64#include <Ifpack2_Parameters.hpp>
65
66namespace Ifpack2 {
67
99template<class GraphType, class KKHandleType>
100class IlukGraph : public Teuchos::Describable {
101public:
102 typedef typename GraphType::local_ordinal_type local_ordinal_type;
103 typedef typename GraphType::global_ordinal_type global_ordinal_type;
104 typedef typename GraphType::node_type node_type;
105
107 typedef Tpetra::RowGraph<local_ordinal_type,
108 global_ordinal_type,
109 node_type> row_graph_type;
111 typedef Tpetra::CrsGraph<local_ordinal_type,
112 global_ordinal_type,
113 node_type> crs_graph_type;
114
115
116
117 typedef typename crs_graph_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
118 typedef typename crs_graph_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
119 typedef typename crs_graph_type::global_inds_host_view_type global_inds_host_view_type;
120 typedef typename crs_graph_type::local_inds_host_view_type local_inds_host_view_type;
121
133 IlukGraph (const Teuchos::RCP<const GraphType>& G,
134 const int levelFill,
135 const int levelOverlap,
136 const double overalloc = 2.);
137
139 virtual ~IlukGraph ();
140
146 void setParameters (const Teuchos::ParameterList& parameterlist);
147
159 void initialize();
160 void initialize(const Teuchos::RCP<KKHandleType>& KernelHandle);
161
163 int getLevelFill () const { return LevelFill_; }
164
166 int getLevelOverlap () const { return LevelOverlap_; }
167
169 Teuchos::RCP<const GraphType> getA_Graph () const {
170 return Graph_;
171 }
172
174 Teuchos::RCP<crs_graph_type> getL_Graph () const {
175 return L_Graph_;
176 }
177
179 Teuchos::RCP<crs_graph_type> getU_Graph () const {
180 return U_Graph_;
181 }
182
184 Teuchos::RCP<const crs_graph_type> getOverlapGraph () const {
185 return OverlapGraph_;
186 }
187
189 size_t getNumGlobalDiagonals() const { return NumGlobalDiagonals_; }
190
191private:
192 typedef typename GraphType::map_type map_type;
193
203
213
215 void constructOverlapGraph();
216
217 Teuchos::RCP<const GraphType> Graph_;
218 Teuchos::RCP<const crs_graph_type> OverlapGraph_;
219 int LevelFill_;
220 int LevelOverlap_;
221 const double Overalloc_;
222 Teuchos::RCP<crs_graph_type> L_Graph_;
223 Teuchos::RCP<crs_graph_type> U_Graph_;
224 size_t NumMyDiagonals_;
225 size_t NumGlobalDiagonals_;
226};
227
228
229template<class GraphType, class KKHandleType>
231IlukGraph (const Teuchos::RCP<const GraphType>& G,
232 const int levelFill,
233 const int levelOverlap,
234 const double overalloc)
235 : Graph_ (G),
236 LevelFill_ (levelFill),
237 LevelOverlap_ (levelOverlap),
238 Overalloc_ (overalloc),
239 NumMyDiagonals_ (0),
240 NumGlobalDiagonals_ (0)
241{
242 TEUCHOS_TEST_FOR_EXCEPTION(Overalloc_ <= 1., std::runtime_error,
243 "Ifpack2::IlukGraph: FATAL: overalloc must be greater than 1.")
244}
245
246
247template<class GraphType, class KKHandleType>
250
251
252template<class GraphType, class KKHandleType>
254setParameters (const Teuchos::ParameterList& parameterlist)
255{
256 getParameter (parameterlist, "iluk level-of-fill", LevelFill_);
257 getParameter (parameterlist, "iluk level-of-overlap", LevelOverlap_);
258}
259
260
261template<class GraphType, class KKHandleType>
262void IlukGraph<GraphType, KKHandleType>::constructOverlapGraph () {
263 // FIXME (mfh 22 Dec 2013) This won't do if we want
264 // RILUK::initialize() to do the right thing (that is,
265 // unconditionally recompute the "symbolic factorization").
266 if (OverlapGraph_ == Teuchos::null) {
267 OverlapGraph_ = createOverlapGraph<GraphType> (Graph_, LevelOverlap_);
268 }
269}
270
271
272template<class GraphType, class KKHandleType>
274{
275 using Teuchos::Array;
276 using Teuchos::ArrayView;
277 using Teuchos::RCP;
278 using Teuchos::rcp;
279 using Teuchos::REDUCE_SUM;
280 using Teuchos::reduceAll;
281
282 size_t NumIn, NumL, NumU;
283 bool DiagFound;
284
285 constructOverlapGraph();
286
287 // Get Maximum Row length
288 const int MaxNumIndices = OverlapGraph_->getLocalMaxNumRowEntries ();
289
290 // FIXME (mfh 23 Dec 2013) Use size_t or whatever
291 // getLocalNumElements() returns, instead of ptrdiff_t.
292 const int NumMyRows = OverlapGraph_->getRowMap ()->getLocalNumElements ();
293
294 using device_type = typename node_type::device_type;
295 using execution_space = typename device_type::execution_space;
296 using dual_view_type = Kokkos::DualView<size_t*,device_type>;
297 dual_view_type numEntPerRow_dv("numEntPerRow",NumMyRows);
298 Tpetra::Details::WrappedDualView<dual_view_type> numEntPerRow(numEntPerRow_dv);
299
300 const auto overalloc = Overalloc_;
301 const auto levelfill = LevelFill_;
302 {
303 // Scoping for the localOverlapGraph access
304 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
305 auto localOverlapGraph = OverlapGraph_->getLocalGraphDevice();
306 Kokkos::parallel_for("CountOverlapGraphRowEntries",
307 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
308 KOKKOS_LAMBDA(const int i)
309 {
310 // Heuristic to get the maximum number of entries per row.
311 int RowMaxNumIndices = localOverlapGraph.rowConst(i).length;
312 numEntPerRow_d(i) = (levelfill == 0) ? RowMaxNumIndices // No additional storage needed
313 : Kokkos::ceil(static_cast<double>(RowMaxNumIndices)
314 * Kokkos::pow(overalloc, levelfill));
315 });
316
317 };
318
319 bool insertError; // No error found yet while inserting entries
320 do {
321 insertError = false;
322 Teuchos::ArrayView<const size_t> a_numEntPerRow(numEntPerRow.getHostView(Tpetra::Access::ReadOnly).data(),NumMyRows);
323 L_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
324 OverlapGraph_->getRowMap (),
325 a_numEntPerRow));
326 U_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
327 OverlapGraph_->getRowMap (),
328 a_numEntPerRow));
329
330 Array<local_ordinal_type> L (MaxNumIndices);
331 Array<local_ordinal_type> U (MaxNumIndices);
332
333 // First we copy the user's graph into L and U, regardless of fill level
334
335 NumMyDiagonals_ = 0;
336
337 for (int i = 0; i< NumMyRows; ++i) {
338 local_inds_host_view_type my_indices;
339 OverlapGraph_->getLocalRowView (i, my_indices);
340
341 // Split into L and U (we don't assume that indices are ordered).
342
343 NumL = 0;
344 NumU = 0;
345 DiagFound = false;
346 NumIn = my_indices.size();
347
348 for (size_t j = 0; j < NumIn; ++j) {
349 const local_ordinal_type k = my_indices[j];
350
351 if (k<NumMyRows) { // Ignore column elements that are not in the square matrix
352
353 if (k==i) {
354 DiagFound = true;
355 }
356 else if (k < i) {
357 L[NumL] = k;
358 NumL++;
359 }
360 else {
361 U[NumU] = k;
362 NumU++;
363 }
364 }
365 }
366
367 // Check in things for this row of L and U
368
369 if (DiagFound) {
370 ++NumMyDiagonals_;
371 }
372 if (NumL) {
373 L_Graph_->insertLocalIndices (i, NumL, L.data());
374 }
375 if (NumU) {
376 U_Graph_->insertLocalIndices (i, NumU, U.data());
377 }
378 }
379
380 if (LevelFill_ > 0) {
381 // Complete Fill steps
382 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
383 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
384 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
385 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
386 RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
387 params->set ("Optimize Storage",false);
388 L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
389 U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
390 L_Graph_->resumeFill ();
391 U_Graph_->resumeFill ();
392
393 // At this point L_Graph and U_Graph are filled with the pattern of input graph,
394 // sorted and have redundant indices (if any) removed. Indices are zero based.
395 // LevelFill is greater than zero, so continue...
396
397 int MaxRC = NumMyRows;
398 std::vector<std::vector<int> > Levels(MaxRC);
399 std::vector<int> LinkList(MaxRC);
400 std::vector<int> CurrentLevel(MaxRC);
401 Array<local_ordinal_type> CurrentRow (MaxRC + 1);
402 std::vector<int> LevelsRowU(MaxRC);
403
404 try {
405 for (int i = 0; i < NumMyRows; ++i) {
406 int First, Next;
407
408 // copy column indices of row into workspace and sort them
409
410 size_t LenL = L_Graph_->getNumEntriesInLocalRow(i);
411 size_t LenU = U_Graph_->getNumEntriesInLocalRow(i);
412 size_t Len = LenL + LenU + 1;
413 CurrentRow.resize(Len);
414 nonconst_local_inds_host_view_type CurrentRow_view(CurrentRow.data(),CurrentRow.size());
415 L_Graph_->getLocalRowCopy(i, CurrentRow_view, LenL); // Get L Indices
416 CurrentRow[LenL] = i; // Put in Diagonal
417 if (LenU > 0) {
418 ArrayView<local_ordinal_type> URowView = CurrentRow.view (LenL+1,LenU);
419 nonconst_local_inds_host_view_type URowView_v(URowView.data(),URowView.size());
420
421 // Get U Indices
422 U_Graph_->getLocalRowCopy (i, URowView_v, LenU);
423 }
424
425 // Construct linked list for current row
426
427 for (size_t j=0; j<Len-1; j++) {
428 LinkList[CurrentRow[j]] = CurrentRow[j+1];
429 CurrentLevel[CurrentRow[j]] = 0;
430 }
431
432 LinkList[CurrentRow[Len-1]] = NumMyRows;
433 CurrentLevel[CurrentRow[Len-1]] = 0;
434
435 // Merge List with rows in U
436
437 First = CurrentRow[0];
438 Next = First;
439 while (Next < i) {
440 int PrevInList = Next;
441 int NextInList = LinkList[Next];
442 int RowU = Next;
443 // Get Indices for this row of U
444 local_inds_host_view_type IndicesU;
445 U_Graph_->getLocalRowView (RowU, IndicesU);
446 // FIXME (mfh 23 Dec 2013) size() returns ptrdiff_t, not int.
447 int LengthRowU = IndicesU.size ();
448
449 int ii;
450
451 // Scan RowU
452
453 for (ii = 0; ii < LengthRowU; /*nop*/) {
454 int CurInList = IndicesU[ii];
455 if (CurInList < NextInList) {
456 // new fill-in
457 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
458 if (NewLevel <= LevelFill_) {
459 LinkList[PrevInList] = CurInList;
460 LinkList[CurInList] = NextInList;
461 PrevInList = CurInList;
462 CurrentLevel[CurInList] = NewLevel;
463 }
464 ii++;
465 }
466 else if (CurInList == NextInList) {
467 PrevInList = NextInList;
468 NextInList = LinkList[PrevInList];
469 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
470 CurrentLevel[CurInList] = std::min (CurrentLevel[CurInList],
471 NewLevel);
472 ii++;
473 }
474 else { // (CurInList > NextInList)
475 PrevInList = NextInList;
476 NextInList = LinkList[PrevInList];
477 }
478 }
479 Next = LinkList[Next];
480 }
481
482 // Put pattern into L and U
483 CurrentRow.resize(0);
484
485 Next = First;
486
487 // Lower
488 while (Next < i) {
489 CurrentRow.push_back(Next);
490 Next = LinkList[Next];
491 }
492
493 // FIXME (mfh 23 Dec 2013) It's not clear to me that
494 // removeLocalIndices works like people expect it to work. In
495 // particular, it does not actually change the column Map.
496 L_Graph_->removeLocalIndices (i); // Delete current set of Indices
497 if (CurrentRow.size() > 0) {
498 L_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
499 }
500
501 // Diagonal
502
503 TEUCHOS_TEST_FOR_EXCEPTION(
504 Next != i, std::runtime_error,
505 "Ifpack2::IlukGraph::initialize: FATAL: U has zero diagonal")
506
507 LevelsRowU[0] = CurrentLevel[Next];
508 Next = LinkList[Next];
509
510 // Upper
511 CurrentRow.resize(0);
512 LenU = 0;
513
514 while (Next < NumMyRows) {
515 LevelsRowU[LenU+1] = CurrentLevel[Next];
516 CurrentRow.push_back (Next);
517 ++LenU;
518 Next = LinkList[Next];
519 }
520
521 // FIXME (mfh 23 Dec 2013) It's not clear to me that
522 // removeLocalIndices works like people expect it to work. In
523 // particular, it does not actually change the column Map.
524
525 U_Graph_->removeLocalIndices (i); // Delete current set of Indices
526 if (LenU > 0) {
527 U_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
528 }
529
530 // Allocate and fill Level info for this row
531 Levels[i] = std::vector<int> (LenU+1);
532 for (size_t jj=0; jj<LenU+1; jj++) {
533 Levels[i][jj] = LevelsRowU[jj];
534 }
535 }
536 }
537 catch (std::runtime_error &e) {
538 insertError = true;
539 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
540 Kokkos::parallel_for("CountOverlapGraphRowEntries",
541 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
542 KOKKOS_LAMBDA(const int i)
543 {
544 const auto numRowEnt = numEntPerRow_d(i);
545 numEntPerRow_d(i) = ceil(static_cast<double>((numRowEnt != 0 ? numRowEnt : 1)) * overalloc);
546 });
547 }
548 const int localInsertError = insertError ? 1 : 0;
549 int globalInsertError = 0;
550 reduceAll (* (OverlapGraph_->getRowMap ()->getComm ()), REDUCE_SUM, 1,
551 &localInsertError, &globalInsertError);
552 insertError = globalInsertError > 0;
553 }
554 } while (insertError); // do until all insertions complete successfully
555
556 // Complete Fill steps
557 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
558 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
559 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
560 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
561 L_Graph_->fillComplete (L_DomainMap, L_RangeMap);//DoOptimizeStorage is default here...
562 U_Graph_->fillComplete (U_DomainMap, U_RangeMap);//DoOptimizeStorage is default here...
563
564 reduceAll<int, size_t> (* (L_DomainMap->getComm ()), REDUCE_SUM, 1,
565 &NumMyDiagonals_, &NumGlobalDiagonals_);
566}
567
568
569template<class GraphType, class KKHandleType>
570void IlukGraph<GraphType, KKHandleType>::initialize(const Teuchos::RCP<KKHandleType>& KernelHandle)
571{
572 using Teuchos::Array;
573 using Teuchos::ArrayView;
574 using Teuchos::RCP;
575 using Teuchos::rcp;
576 using Teuchos::REDUCE_SUM;
577 using Teuchos::reduceAll;
578
579 typedef typename crs_graph_type::local_graph_device_type local_graph_device_type;
580 typedef typename local_graph_device_type::size_type size_type;
581 typedef typename local_graph_device_type::data_type data_type;
582 typedef typename local_graph_device_type::array_layout array_layout;
583 typedef typename local_graph_device_type::device_type device_type;
584
585 typedef typename Kokkos::View<size_type*, array_layout, device_type> lno_row_view_t;
586 typedef typename Kokkos::View<data_type*, array_layout, device_type> lno_nonzero_view_t;
587
588 constructOverlapGraph();
589
590 // FIXME (mfh 23 Dec 2013) Use size_t or whatever
591 // getLocalNumElements() returns, instead of ptrdiff_t.
592 const int NumMyRows = OverlapGraph_->getRowMap()->getLocalNumElements();
593 auto localOverlapGraph = OverlapGraph_->getLocalGraphDevice();
594
595 if (KernelHandle->get_spiluk_handle()->get_nrows() < static_cast<size_type>(NumMyRows)) {
596 KernelHandle->get_spiluk_handle()->reset_handle(NumMyRows,
597 KernelHandle->get_spiluk_handle()->get_nnzL(),
598 KernelHandle->get_spiluk_handle()->get_nnzU());
599 }
600
601 lno_row_view_t L_row_map("L_row_map", NumMyRows + 1);
602 lno_nonzero_view_t L_entries("L_entries", KernelHandle->get_spiluk_handle()->get_nnzL());
603 lno_row_view_t U_row_map("U_row_map", NumMyRows + 1);
604 lno_nonzero_view_t U_entries("U_entries", KernelHandle->get_spiluk_handle()->get_nnzU());
605
606 bool symbolicError;
607 do {
608 symbolicError = false;
609 try {
610 KokkosSparse::Experimental::spiluk_symbolic( KernelHandle.getRawPtr(), LevelFill_,
611 localOverlapGraph.row_map, localOverlapGraph.entries,
612 L_row_map, L_entries, U_row_map, U_entries );
613 }
614 catch (std::runtime_error &e) {
615 symbolicError = true;
616 data_type nnzL = static_cast<data_type>(Overalloc_)*L_entries.extent(0);
617 data_type nnzU = static_cast<data_type>(Overalloc_)*U_entries.extent(0);
618 KernelHandle->get_spiluk_handle()->reset_handle(NumMyRows, nnzL, nnzU);
619 Kokkos::resize(L_entries, KernelHandle->get_spiluk_handle()->get_nnzL());
620 Kokkos::resize(U_entries, KernelHandle->get_spiluk_handle()->get_nnzU());
621 }
622 const int localSymbolicError = symbolicError ? 1 : 0;
623 int globalSymbolicError = 0;
624 reduceAll (* (OverlapGraph_->getRowMap ()->getComm ()), REDUCE_SUM, 1,
625 &localSymbolicError, &globalSymbolicError);
626 symbolicError = globalSymbolicError > 0;
627 } while (symbolicError);
628
629 Kokkos::resize(L_entries, KernelHandle->get_spiluk_handle()->get_nnzL());
630 Kokkos::resize(U_entries, KernelHandle->get_spiluk_handle()->get_nnzU());
631
632 RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
633 params->set ("Optimize Storage",false);
634
635 L_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
636 OverlapGraph_->getRowMap (),
637 L_row_map, L_entries));
638 U_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
639 OverlapGraph_->getRowMap (),
640 U_row_map, U_entries));
641
642 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
643 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
644 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
645 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
646 L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
647 U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
648}
649
650} // namespace Ifpack2
651
652#endif /* IFPACK2_ILUK_GRAPH_HPP */
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition Ifpack2_IlukGraph.hpp:100
size_t getNumGlobalDiagonals() const
Returns the global number of diagonals in the ILU(k) graph.
Definition Ifpack2_IlukGraph.hpp:189
Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > crs_graph_type
Tpetra::CrsGraph specialization used by this class.
Definition Ifpack2_IlukGraph.hpp:113
Teuchos::RCP< crs_graph_type > getU_Graph() const
Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph.
Definition Ifpack2_IlukGraph.hpp:179
int getLevelFill() const
The level of fill used to construct this graph.
Definition Ifpack2_IlukGraph.hpp:163
int getLevelOverlap() const
The level of overlap used to construct this graph.
Definition Ifpack2_IlukGraph.hpp:166
Teuchos::RCP< const GraphType > getA_Graph() const
Returns the original graph given.
Definition Ifpack2_IlukGraph.hpp:169
void initialize()
Set up the graph structure of the L and U factors.
Definition Ifpack2_IlukGraph.hpp:273
Teuchos::RCP< const crs_graph_type > getOverlapGraph() const
Returns the the overlapped graph.
Definition Ifpack2_IlukGraph.hpp:184
void setParameters(const Teuchos::ParameterList &parameterlist)
Set parameters.
Definition Ifpack2_IlukGraph.hpp:254
Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > row_graph_type
Tpetra::RowGraph specialization used by this class.
Definition Ifpack2_IlukGraph.hpp:109
Teuchos::RCP< crs_graph_type > getL_Graph() const
Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph.
Definition Ifpack2_IlukGraph.hpp:174
IlukGraph(const Teuchos::RCP< const GraphType > &G, const int levelFill, const int levelOverlap, const double overalloc=2.)
Constructor.
Definition Ifpack2_IlukGraph.hpp:231
virtual ~IlukGraph()
IlukGraph Destructor.
Definition Ifpack2_IlukGraph.hpp:248
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74
void getParameter(const Teuchos::ParameterList &params, const std::string &name, T &value)
Set a value from a ParameterList if a parameter with the specified name exists.
Definition Ifpack2_Parameters.hpp:59
Teuchos::RCP< const GraphType > createOverlapGraph(const Teuchos::RCP< const GraphType > &inputGraph, const int overlapLevel)
Construct an overlapped graph for use with Ifpack2 preconditioners.
Definition Ifpack2_CreateOverlapGraph.hpp:73