Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_RILUK_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_CRSRILUK_DEF_HPP
11#define IFPACK2_CRSRILUK_DEF_HPP
12
14#include "Ifpack2_LocalFilter.hpp"
15#include "Tpetra_CrsMatrix.hpp"
16#include "Teuchos_StandardParameterEntryValidators.hpp"
17#include "Ifpack2_LocalSparseTriangularSolver.hpp"
18#include "Ifpack2_Details_getParamTryingTypes.hpp"
19#include "Ifpack2_Details_getCrsMatrix.hpp"
20#include "Kokkos_Sort.hpp"
21#include "KokkosSparse_spiluk.hpp"
22#include "KokkosSparse_Utils.hpp"
23#include "KokkosKernels_Sorting.hpp"
24#include "KokkosSparse_IOUtils.hpp"
25
26namespace Ifpack2 {
27
28namespace Details {
29struct IlukImplType {
30 enum Enum {
31 Serial,
32 KSPILUK
33 };
34
35 static void loadPLTypeOption(Teuchos::Array<std::string>& type_strs, Teuchos::Array<Enum>& type_enums) {
36 type_strs.resize(2);
37 type_strs[0] = "Serial";
38 type_strs[1] = "KSPILUK";
39 type_enums.resize(2);
40 type_enums[0] = Serial;
41 type_enums[1] = KSPILUK;
42 }
43};
44} // namespace Details
45
46template <class MatrixType>
47RILUK<MatrixType>::RILUK(const Teuchos::RCP<const row_matrix_type>& Matrix_in, const Teuchos::RCP<const coord_type>& Matrix_in_coordinates)
48 : A_(Matrix_in)
49 , A_coordinates_(Matrix_in_coordinates)
50 , LevelOfFill_(0)
51 , Overalloc_(2.)
52 , isAllocated_(false)
53 , isInitialized_(false)
54 , isComputed_(false)
55 , numInitialize_(0)
56 , numCompute_(0)
57 , numApply_(0)
58 , initializeTime_(0.0)
59 , computeTime_(0.0)
60 , applyTime_(0.0)
61 , RelaxValue_(Teuchos::ScalarTraits<magnitude_type>::zero())
62 , Athresh_(Teuchos::ScalarTraits<magnitude_type>::zero())
63 , Rthresh_(Teuchos::ScalarTraits<magnitude_type>::one())
64 , isKokkosKernelsSpiluk_(false)
65 , isKokkosKernelsStream_(false)
66 , num_streams_(0)
67 , hasStreamReordered_(false)
68 , hasStreamsWithRCB_(false) {
69 allocateSolvers();
70}
71
72template <class MatrixType>
73RILUK<MatrixType>::RILUK(const Teuchos::RCP<const crs_matrix_type>& Matrix_in, const Teuchos::RCP<const coord_type>& Matrix_in_coordinates)
74 : A_(Matrix_in)
75 , A_coordinates_(Matrix_in_coordinates)
76 , LevelOfFill_(0)
77 , Overalloc_(2.)
78 , isAllocated_(false)
79 , isInitialized_(false)
80 , isComputed_(false)
81 , numInitialize_(0)
82 , numCompute_(0)
83 , numApply_(0)
84 , initializeTime_(0.0)
85 , computeTime_(0.0)
86 , applyTime_(0.0)
87 , RelaxValue_(Teuchos::ScalarTraits<magnitude_type>::zero())
88 , Athresh_(Teuchos::ScalarTraits<magnitude_type>::zero())
89 , Rthresh_(Teuchos::ScalarTraits<magnitude_type>::one())
90 , isKokkosKernelsSpiluk_(false)
91 , isKokkosKernelsStream_(false)
92 , num_streams_(0)
93 , hasStreamReordered_(false)
94 , hasStreamsWithRCB_(false) {
95 allocateSolvers();
96}
97
98template <class MatrixType>
100 if (!isKokkosKernelsStream_) {
101 if (Teuchos::nonnull(KernelHandle_)) {
102 KernelHandle_->destroy_spiluk_handle();
103 }
104 } else {
105 for (size_t i = 0; i < KernelHandle_v_.size(); i++) {
106 if (Teuchos::nonnull(KernelHandle_v_[i])) {
107 KernelHandle_v_[i]->destroy_spiluk_handle();
108 }
109 }
110 }
111}
112
113template <class MatrixType>
115 L_solver_ = Teuchos::rcp(new LocalSparseTriangularSolver<row_matrix_type>());
116 L_solver_->setObjectLabel("lower");
117 U_solver_ = Teuchos::rcp(new LocalSparseTriangularSolver<row_matrix_type>());
118 U_solver_->setObjectLabel("upper");
119}
120
121template <class MatrixType>
122void RILUK<MatrixType>::setMatrix(const Teuchos::RCP<const row_matrix_type>& A) {
123 // It's legal for A to be null; in that case, you may not call
124 // initialize() until calling setMatrix() with a nonnull input.
125 // Regardless, setting the matrix invalidates any previous
126 // factorization.
127 if (A.getRawPtr() != A_.getRawPtr()) {
128 isAllocated_ = false;
129 isInitialized_ = false;
130 isComputed_ = false;
131 A_local_ = Teuchos::null;
132 Graph_ = Teuchos::null;
133
134 // The sparse triangular solvers get a triangular factor as their
135 // input matrix. The triangular factors L_ and U_ are getting
136 // reset, so we reset the solvers' matrices to null. Do that
137 // before setting L_ and U_ to null, so that latter step actually
138 // frees the factors.
139 if (!L_solver_.is_null()) {
140 L_solver_->setMatrix(Teuchos::null);
141 }
142 if (!U_solver_.is_null()) {
143 U_solver_->setMatrix(Teuchos::null);
144 }
145
146 L_ = Teuchos::null;
147 U_ = Teuchos::null;
148 D_ = Teuchos::null;
149 A_ = A;
150 }
151}
152
153template <class MatrixType>
154void RILUK<MatrixType>::setCoord(const Teuchos::RCP<const coord_type>& A_coordinates) {
155 if (A_coordinates.getRawPtr() != A_coordinates_.getRawPtr()) {
156 A_coordinates_ = A_coordinates;
157 }
158}
159
160template <class MatrixType>
163 TEUCHOS_TEST_FOR_EXCEPTION(
164 L_.is_null(), std::runtime_error,
165 "Ifpack2::RILUK::getL: The L factor "
166 "is null. Please call initialize() (and preferably also compute()) "
167 "before calling this method. If the input matrix has not yet been set, "
168 "you must first call setMatrix() with a nonnull input matrix before you "
169 "may call initialize() or compute().");
170 return *L_;
171}
172
173template <class MatrixType>
174const Tpetra::Vector<typename RILUK<MatrixType>::scalar_type,
179 TEUCHOS_TEST_FOR_EXCEPTION(
180 D_.is_null(), std::runtime_error,
181 "Ifpack2::RILUK::getD: The D factor "
182 "(of diagonal entries) is null. Please call initialize() (and "
183 "preferably also compute()) before calling this method. If the input "
184 "matrix has not yet been set, you must first call setMatrix() with a "
185 "nonnull input matrix before you may call initialize() or compute().");
186 return *D_;
187}
188
189template <class MatrixType>
192 TEUCHOS_TEST_FOR_EXCEPTION(
193 U_.is_null(), std::runtime_error,
194 "Ifpack2::RILUK::getU: The U factor "
195 "is null. Please call initialize() (and preferably also compute()) "
196 "before calling this method. If the input matrix has not yet been set, "
197 "you must first call setMatrix() with a nonnull input matrix before you "
198 "may call initialize() or compute().");
199 return *U_;
200}
201
202template <class MatrixType>
204 TEUCHOS_TEST_FOR_EXCEPTION(
205 A_.is_null(), std::runtime_error,
206 "Ifpack2::RILUK::getNodeSmootherComplexity: "
207 "The input matrix A is null. Please call setMatrix() with a nonnull "
208 "input matrix, then call compute(), before calling this method.");
209 // RILUK methods cost roughly one apply + the nnz in the upper+lower triangles
210 if (!L_.is_null() && !U_.is_null())
211 return A_->getLocalNumEntries() + L_->getLocalNumEntries() + U_->getLocalNumEntries();
212 else
213 return 0;
214}
215
216template <class MatrixType>
217Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
221 TEUCHOS_TEST_FOR_EXCEPTION(
222 A_.is_null(), std::runtime_error,
223 "Ifpack2::RILUK::getDomainMap: "
224 "The matrix is null. Please call setMatrix() with a nonnull input "
225 "before calling this method.");
226
227 // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
228 TEUCHOS_TEST_FOR_EXCEPTION(
229 Graph_.is_null(), std::runtime_error,
230 "Ifpack2::RILUK::getDomainMap: "
231 "The computed graph is null. Please call initialize() before calling "
232 "this method.");
233 return Graph_->getL_Graph()->getDomainMap();
234}
235
236template <class MatrixType>
237Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
241 TEUCHOS_TEST_FOR_EXCEPTION(
242 A_.is_null(), std::runtime_error,
243 "Ifpack2::RILUK::getRangeMap: "
244 "The matrix is null. Please call setMatrix() with a nonnull input "
245 "before calling this method.");
246
247 // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
248 TEUCHOS_TEST_FOR_EXCEPTION(
249 Graph_.is_null(), std::runtime_error,
250 "Ifpack2::RILUK::getRangeMap: "
251 "The computed graph is null. Please call initialize() before calling "
252 "this method.");
253 return Graph_->getL_Graph()->getRangeMap();
254}
255
256template <class MatrixType>
258 using Teuchos::null;
259 using Teuchos::rcp;
260
261 if (!isAllocated_) {
262 if (!isKokkosKernelsStream_) {
263 // Deallocate any existing storage. This avoids storing 2x
264 // memory, since RCP op= does not deallocate until after the
265 // assignment.
266 L_ = null;
267 U_ = null;
268 D_ = null;
269
270 // Allocate Matrix using ILUK graphs
271 L_ = rcp(new crs_matrix_type(Graph_->getL_Graph()));
272 U_ = rcp(new crs_matrix_type(Graph_->getU_Graph()));
273 L_->setAllToScalar(STS::zero()); // Zero out L and U matrices
274 U_->setAllToScalar(STS::zero());
275
276 // FIXME (mfh 24 Jan 2014) This assumes domain == range Map for L and U.
277 L_->fillComplete();
278 U_->fillComplete();
279 D_ = rcp(new vec_type(Graph_->getL_Graph()->getRowMap()));
280 } else {
281 L_v_ = std::vector<Teuchos::RCP<crs_matrix_type> >(num_streams_, null);
282 U_v_ = std::vector<Teuchos::RCP<crs_matrix_type> >(num_streams_, null);
283 for (int i = 0; i < num_streams_; i++) {
284 L_v_[i] = rcp(new crs_matrix_type(Graph_v_[i]->getL_Graph()));
285 U_v_[i] = rcp(new crs_matrix_type(Graph_v_[i]->getU_Graph()));
286 L_v_[i]->setAllToScalar(STS::zero()); // Zero out L and U matrices
287 U_v_[i]->setAllToScalar(STS::zero());
288
289 L_v_[i]->fillComplete();
290 U_v_[i]->fillComplete();
291 }
292 }
293 }
294 isAllocated_ = true;
295}
296
297template <class MatrixType>
299 setParameters(const Teuchos::ParameterList& params) {
300 using Details::getParamTryingTypes;
301 using Teuchos::Array;
302 using Teuchos::ParameterList;
303 using Teuchos::RCP;
304 const char prefix[] = "Ifpack2::RILUK: ";
305
306 // Default values of the various parameters.
307 int fillLevel = 0;
308 magnitude_type absThresh = STM::zero();
309 magnitude_type relThresh = STM::one();
310 magnitude_type relaxValue = STM::zero();
311 double overalloc = 2.;
312 int nstreams = 0;
313
314 // "fact: iluk level-of-fill" parsing is more complicated, because
315 // we want to allow as many types as make sense. int is the native
316 // type, but we also want to accept double (for backwards
317 // compatibilty with ILUT). You can't cast arbitrary magnitude_type
318 // (e.g., Sacado::MP::Vector) to int, so we use float instead, to
319 // get coverage of the most common magnitude_type cases. Weirdly,
320 // there's an Ifpack2 test that sets the fill level as a
321 // global_ordinal_type.
322 {
323 const std::string paramName("fact: iluk level-of-fill");
324 getParamTryingTypes<int, int, global_ordinal_type, double, float>(fillLevel, params, paramName, prefix);
326 // For the other parameters, we prefer magnitude_type, but allow
327 // double for backwards compatibility.
329 const std::string paramName("fact: absolute threshold");
330 getParamTryingTypes<magnitude_type, magnitude_type, double>(absThresh, params, paramName, prefix);
331 }
332 {
333 const std::string paramName("fact: relative threshold");
334 getParamTryingTypes<magnitude_type, magnitude_type, double>(relThresh, params, paramName, prefix);
335 }
336 {
337 const std::string paramName("fact: relax value");
338 getParamTryingTypes<magnitude_type, magnitude_type, double>(relaxValue, params, paramName, prefix);
339 }
340 {
341 const std::string paramName("fact: iluk overalloc");
342 getParamTryingTypes<double, double>(overalloc, params, paramName, prefix);
343 }
344
345 // Parsing implementation type
346 Details::IlukImplType::Enum ilukimplType = Details::IlukImplType::Serial;
347 do {
348 static const char typeName[] = "fact: type";
349
350 if (!params.isType<std::string>(typeName)) break;
351
352 // Map std::string <-> IlukImplType::Enum.
353 Array<std::string> ilukimplTypeStrs;
354 Array<Details::IlukImplType::Enum> ilukimplTypeEnums;
355 Details::IlukImplType::loadPLTypeOption(ilukimplTypeStrs, ilukimplTypeEnums);
356 Teuchos::StringToIntegralParameterEntryValidator<Details::IlukImplType::Enum>
357 s2i(ilukimplTypeStrs(), ilukimplTypeEnums(), typeName, false);
358
359 ilukimplType = s2i.getIntegralValue(params.get<std::string>(typeName));
360 } while (0);
361
362 if (ilukimplType == Details::IlukImplType::KSPILUK) {
363 this->isKokkosKernelsSpiluk_ = true;
364 } else {
365 this->isKokkosKernelsSpiluk_ = false;
366 }
367
368 {
369 const std::string paramName("fact: kspiluk number-of-streams");
370 getParamTryingTypes<int, int, global_ordinal_type>(nstreams, params, paramName, prefix);
371 }
372
373 // Forward to trisolvers.
374 L_solver_->setParameters(params);
375 U_solver_->setParameters(params);
377 // "Commit" the values only after validating all of them. This
378 // ensures that there are no side effects if this routine throws an
379 // exception.
380
381 LevelOfFill_ = fillLevel;
382 Overalloc_ = overalloc;
383#ifdef KOKKOS_ENABLE_OPENMP
384 if constexpr (std::is_same_v<execution_space, Kokkos::OpenMP>) {
385 nstreams = std::min(nstreams, execution_space{}.concurrency());
386 }
387#endif
388 num_streams_ = nstreams;
389
390 if (num_streams_ >= 1) {
391 this->isKokkosKernelsStream_ = true;
392 // Will we do reordering in streams?
393 if (params.isParameter("fact: kspiluk reordering in streams"))
394 hasStreamReordered_ = params.get<bool>("fact: kspiluk reordering in streams");
395
396 if (params.isParameter("fact: kspiluk streams use rcb"))
397 hasStreamsWithRCB_ = params.get<bool>("fact: kspiluk streams use rcb");
398 } else {
399 this->isKokkosKernelsStream_ = false;
400 }
401
402 // mfh 28 Nov 2012: The previous code would not assign Athresh_,
403 // Rthresh_, or RelaxValue_, if the read-in value was -1. I don't
404 // know if keeping this behavior is correct, but I'll keep it just
405 // so as not to change previous behavior.
406
407 if (absThresh != -STM::one()) {
408 Athresh_ = absThresh;
409 }
410 if (relThresh != -STM::one()) {
411 Rthresh_ = relThresh;
412 }
413 if (relaxValue != -STM::one()) {
414 RelaxValue_ = relaxValue;
415 }
416}
418template <class MatrixType>
419Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
421 return Teuchos::rcp_implicit_cast<const row_matrix_type>(A_);
422}
423
424template <class MatrixType>
425Teuchos::RCP<const typename RILUK<MatrixType>::crs_matrix_type>
427 return Teuchos::rcp_dynamic_cast<const crs_matrix_type>(A_, true);
428}
430template <class MatrixType>
431Teuchos::RCP<const typename RILUK<MatrixType>::coord_type>
433 return A_coordinates_;
434}
435
436template <class MatrixType>
437Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
438RILUK<MatrixType>::makeLocalFilter(const Teuchos::RCP<const row_matrix_type>& A) {
439 using Teuchos::RCP;
440 using Teuchos::rcp;
441 using Teuchos::rcp_dynamic_cast;
442 using Teuchos::rcp_implicit_cast;
443
444 // If A_'s communicator only has one process, or if its column and
445 // row Maps are the same, then it is already local, so use it
446 // directly.
447 if (A->getRowMap()->getComm()->getSize() == 1 ||
448 A->getRowMap()->isSameAs(*(A->getColMap()))) {
449 return A;
450 }
451
452 // If A_ is already a LocalFilter, then use it directly. This
453 // should be the case if RILUK is being used through
454 // AdditiveSchwarz, for example.
455 RCP<const LocalFilter<row_matrix_type> > A_lf_r =
456 rcp_dynamic_cast<const LocalFilter<row_matrix_type> >(A);
457 if (!A_lf_r.is_null()) {
458 return rcp_implicit_cast<const row_matrix_type>(A_lf_r);
459 } else {
460 // A_'s communicator has more than one process, its row Map and
461 // its column Map differ, and A_ is not a LocalFilter. Thus, we
462 // have to wrap it in a LocalFilter.
463 return rcp(new LocalFilter<row_matrix_type>(A));
464 }
465}
466
467template <class MatrixType>
469 using Teuchos::Array;
470 using Teuchos::ArrayView;
471 using Teuchos::RCP;
472 using Teuchos::rcp;
473 using Teuchos::rcp_const_cast;
474 using Teuchos::rcp_dynamic_cast;
475 using Teuchos::rcp_implicit_cast;
476 typedef Tpetra::CrsGraph<local_ordinal_type,
478 node_type>
479 crs_graph_type;
480 typedef Tpetra::Map<local_ordinal_type,
482 node_type>
483 crs_map_type;
484 const char prefix[] = "Ifpack2::RILUK::initialize: ";
485
486 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix is null. Please "
487 "call setMatrix() with a nonnull input before calling this method.");
488 TEUCHOS_TEST_FOR_EXCEPTION(!A_->isFillComplete(), std::runtime_error, prefix << "The matrix is not "
489 "fill complete. You may not invoke initialize() or compute() with this "
490 "matrix until the matrix is fill complete. If your matrix is a "
491 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
492 "range Maps, if appropriate) before calling this method.");
493
494 Teuchos::Time timer("RILUK::initialize");
495 double startTime = timer.wallTime();
496 { // Start timing
497 Teuchos::TimeMonitor timeMon(timer);
498
499 // Calling initialize() means that the user asserts that the graph
500 // of the sparse matrix may have changed. We must not just reuse
501 // the previous graph in that case.
502 //
503 // Regarding setting isAllocated_ to false: Eventually, we may want
504 // some kind of clever memory reuse strategy, but it's always
505 // correct just to blow everything away and start over.
506 isInitialized_ = false;
507 isAllocated_ = false;
508 isComputed_ = false;
509 Graph_ = Teuchos::null;
510 Y_tmp_ = nullptr;
511 reordered_x_ = nullptr;
512 reordered_y_ = nullptr;
513
514 if (isKokkosKernelsStream_) {
515 Graph_v_ = std::vector<Teuchos::RCP<iluk_graph_type> >(num_streams_);
516 A_local_diagblks_v_ = std::vector<local_matrix_device_type>(num_streams_);
517 for (int i = 0; i < num_streams_; i++) {
518 Graph_v_[i] = Teuchos::null;
519 }
520 }
521
522 A_local_ = makeLocalFilter(A_);
523
524 TEUCHOS_TEST_FOR_EXCEPTION(
525 A_local_.is_null(), std::logic_error,
526 "Ifpack2::RILUK::initialize: "
527 "makeLocalFilter returned null; it failed to compute A_local. "
528 "Please report this bug to the Ifpack2 developers.");
529
530 // FIXME (mfh 24 Jan 2014, 26 Mar 2014) It would be more efficient
531 // to rewrite RILUK so that it works with any RowMatrix input, not
532 // just CrsMatrix. (That would require rewriting IlukGraph to
533 // handle a Tpetra::RowGraph.) However, to make it work for now,
534 // we just copy the input matrix if it's not a CrsMatrix.
535
536 {
537 A_local_crs_ = Details::getCrsMatrix(A_local_);
538 if (A_local_crs_.is_null()) {
539 local_ordinal_type numRows = A_local_->getLocalNumRows();
540 Array<size_t> entriesPerRow(numRows);
541 for (local_ordinal_type i = 0; i < numRows; i++) {
542 entriesPerRow[i] = A_local_->getNumEntriesInLocalRow(i);
543 }
544 A_local_crs_nc_ =
545 rcp(new crs_matrix_type(A_local_->getRowMap(),
546 A_local_->getColMap(),
547 entriesPerRow()));
548 // copy entries into A_local_crs
549 nonconst_local_inds_host_view_type indices("indices", A_local_->getLocalMaxNumRowEntries());
550 nonconst_values_host_view_type values("values", A_local_->getLocalMaxNumRowEntries());
551 for (local_ordinal_type i = 0; i < numRows; i++) {
552 size_t numEntries = 0;
553 A_local_->getLocalRowCopy(i, indices, values, numEntries);
554 A_local_crs_nc_->insertLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()), indices.data());
555 }
556 A_local_crs_nc_->fillComplete(A_local_->getDomainMap(), A_local_->getRangeMap());
557 A_local_crs_ = rcp_const_cast<const crs_matrix_type>(A_local_crs_nc_);
559 if (!isKokkosKernelsStream_) {
560 Graph_ = rcp(new Ifpack2::IlukGraph<crs_graph_type, kk_handle_type>(A_local_crs_->getCrsGraph(),
561 LevelOfFill_, 0, Overalloc_));
562 } else {
563 std::vector<int> weights(num_streams_);
564 std::fill(weights.begin(), weights.end(), 1);
565 exec_space_instances_ = Kokkos::Experimental::partition_space(execution_space(), weights);
566
567 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
568 if (!hasStreamReordered_) {
569 if (hasStreamsWithRCB_) {
570 TEUCHOS_TEST_FOR_EXCEPTION(A_coordinates_.is_null(), std::runtime_error, prefix << "The coordinates associated with rows of the input matrix is null while RILUK uses streams with RCB. Please call setCoord() with a nonnull input before calling this method.");
571 auto A_coordinates_lcl = A_coordinates_->getLocalViewDevice(Tpetra::Access::ReadOnly);
572 perm_rcb_ = perm_view_t(Kokkos::view_alloc(Kokkos::WithoutInitializing, "perm_rcb_"), A_coordinates_lcl.extent(0));
573 coors_rcb_ = coors_view_t(Kokkos::view_alloc(Kokkos::WithoutInitializing, "coors_rcb_"), A_coordinates_lcl.extent(0), A_coordinates_lcl.extent(1));
574 Kokkos::deep_copy(coors_rcb_, A_coordinates_lcl);
575 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, coors_rcb_,
576 A_local_diagblks_v_, perm_rcb_);
577 } else {
578 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_);
579 }
580 } else {
581 perm_v_ = KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, true);
582 reverse_perm_v_.resize(perm_v_.size());
583 for (size_t istream = 0; istream < perm_v_.size(); ++istream) {
584 using perm_type = typename lno_nonzero_view_t::non_const_type;
585 const auto perm = perm_v_[istream];
586 const auto perm_length = perm.extent(0);
587 perm_type reverse_perm(
588 Kokkos::view_alloc(Kokkos::WithoutInitializing, "reverse_perm"),
589 perm_length);
590 Kokkos::parallel_for(
591 Kokkos::RangePolicy<execution_space>(exec_space_instances_[istream], 0, perm_length),
592 KOKKOS_LAMBDA(const local_ordinal_type ii) {
593 reverse_perm(perm(ii)) = ii;
594 });
595 reverse_perm_v_[istream] = reverse_perm;
596 }
597 }
598
599 A_local_diagblks_rowmap_v_ = std::vector<lno_row_view_t>(num_streams_);
600 A_local_diagblks_entries_v_ = std::vector<lno_nonzero_view_t>(num_streams_);
601 A_local_diagblks_values_v_ = std::vector<scalar_nonzero_view_t>(num_streams_);
602
603 for (int i = 0; i < num_streams_; i++) {
604 A_local_diagblks_rowmap_v_[i] = A_local_diagblks_v_[i].graph.row_map;
605 A_local_diagblks_entries_v_[i] = A_local_diagblks_v_[i].graph.entries;
606 A_local_diagblks_values_v_[i] = A_local_diagblks_v_[i].values;
607
608 Teuchos::RCP<const crs_map_type> A_local_diagblks_RowMap = rcp(new crs_map_type(A_local_diagblks_v_[i].numRows(),
609 A_local_diagblks_v_[i].numRows(),
610 A_local_crs_->getRowMap()->getComm()));
611 Teuchos::RCP<const crs_map_type> A_local_diagblks_ColMap = rcp(new crs_map_type(A_local_diagblks_v_[i].numCols(),
612 A_local_diagblks_v_[i].numCols(),
613 A_local_crs_->getColMap()->getComm()));
614 Teuchos::RCP<crs_matrix_type> A_local_diagblks = rcp(new crs_matrix_type(A_local_diagblks_RowMap,
615 A_local_diagblks_ColMap,
616 A_local_diagblks_v_[i]));
617 Graph_v_[i] = rcp(new Ifpack2::IlukGraph<crs_graph_type, kk_handle_type>(A_local_diagblks->getCrsGraph(),
618 LevelOfFill_, 0, Overalloc_));
619 }
620 }
621 }
622
623 if (this->isKokkosKernelsSpiluk_) {
624 if (!isKokkosKernelsStream_) {
625 this->KernelHandle_ = Teuchos::rcp(new kk_handle_type());
626 KernelHandle_->create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1,
627 A_local_->getLocalNumRows(),
628 2 * A_local_->getLocalNumEntries() * (LevelOfFill_ + 1),
629 2 * A_local_->getLocalNumEntries() * (LevelOfFill_ + 1));
630 Graph_->initialize(KernelHandle_); // this calls spiluk_symbolic
631 } else {
632 KernelHandle_v_ = std::vector<Teuchos::RCP<kk_handle_type> >(num_streams_);
633 for (int i = 0; i < num_streams_; i++) {
634 KernelHandle_v_[i] = Teuchos::rcp(new kk_handle_type());
635 KernelHandle_v_[i]->create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1,
636 A_local_diagblks_v_[i].numRows(),
637 2 * A_local_diagblks_v_[i].nnz() * (LevelOfFill_ + 1),
638 2 * A_local_diagblks_v_[i].nnz() * (LevelOfFill_ + 1));
639 Graph_v_[i]->initialize(KernelHandle_v_[i]); // this calls spiluk_symbolic
640 }
641 }
642 } else {
643 Graph_->initialize();
644 }
645
646 allocate_L_and_U();
647 checkOrderingConsistency(*A_local_);
648 if (!isKokkosKernelsStream_) {
649 L_solver_->setMatrix(L_);
650 } else {
651 L_solver_->setStreamInfo(isKokkosKernelsStream_, num_streams_, exec_space_instances_);
652 L_solver_->setMatrices(L_v_);
653 }
654 L_solver_->initialize();
655
656 if (!isKokkosKernelsStream_) {
657 U_solver_->setMatrix(U_);
658 } else {
659 U_solver_->setStreamInfo(isKokkosKernelsStream_, num_streams_, exec_space_instances_);
660 U_solver_->setMatrices(U_v_);
661 }
662 U_solver_->initialize();
663
664 // Do not call initAllValues. compute() always calls initAllValues to
665 // fill L and U with possibly new numbers. initialize() is concerned
666 // only with the nonzero pattern.
667 } // Stop timing
668
669 isInitialized_ = true;
670 ++numInitialize_;
671 initializeTime_ += (timer.wallTime() - startTime);
672}
673
674template <class MatrixType>
675void RILUK<MatrixType>::
676 checkOrderingConsistency(const row_matrix_type& A) {
677 // First check that the local row map ordering is the same as the local portion of the column map.
678 // The extraction of the strictly lower/upper parts of A, as well as the factorization,
679 // implicitly assume that this is the case.
680 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getLocalElementList();
681 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getLocalElementList();
682 bool gidsAreConsistentlyOrdered = true;
683 global_ordinal_type indexOfInconsistentGID = 0;
684 for (global_ordinal_type i = 0; i < rowGIDs.size(); ++i) {
685 if (rowGIDs[i] != colGIDs[i]) {
686 gidsAreConsistentlyOrdered = false;
687 indexOfInconsistentGID = i;
688 break;
689 }
690 }
691 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered == false, std::runtime_error,
692 "The ordering of the local GIDs in the row and column maps is not the same"
693 << std::endl
694 << "at index " << indexOfInconsistentGID
695 << ". Consistency is required, as all calculations are done with"
696 << std::endl
697 << "local indexing.");
698}
699
700template <class MatrixType>
701void RILUK<MatrixType>::
702 initAllValues(const row_matrix_type& A) {
703 using Teuchos::ArrayRCP;
704 using Teuchos::Comm;
705 using Teuchos::ptr;
706 using Teuchos::RCP;
707 using Teuchos::REDUCE_SUM;
708 using Teuchos::reduceAll;
709 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
710
711 size_t NumIn = 0, NumL = 0, NumU = 0;
712 bool DiagFound = false;
713 size_t NumNonzeroDiags = 0;
714 size_t MaxNumEntries = A.getGlobalMaxNumRowEntries();
715
716 // Allocate temporary space for extracting the strictly
717 // lower and upper parts of the matrix A.
718 nonconst_local_inds_host_view_type InI("InI", MaxNumEntries);
719 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
720 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
721 nonconst_values_host_view_type InV("InV", MaxNumEntries);
722 Teuchos::Array<scalar_type> LV(MaxNumEntries);
723 Teuchos::Array<scalar_type> UV(MaxNumEntries);
724
725 // Check if values should be inserted or replaced
726 const bool ReplaceValues = L_->isStaticGraph() || L_->isLocallyIndexed();
727
728 L_->resumeFill();
729 U_->resumeFill();
730 if (ReplaceValues) {
731 L_->setAllToScalar(STS::zero()); // Zero out L and U matrices
732 U_->setAllToScalar(STS::zero());
733 }
734
735 D_->putScalar(STS::zero()); // Set diagonal values to zero
736 auto DV = Kokkos::subview(D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
737
738 RCP<const map_type> rowMap = L_->getRowMap();
739
740 // First we copy the user's matrix into L and U, regardless of fill level.
741 // It is important to note that L and U are populated using local indices.
742 // This means that if the row map GIDs are not monotonically increasing
743 // (i.e., permuted or gappy), then the strictly lower (upper) part of the
744 // matrix is not the one that you would get if you based L (U) on GIDs.
745 // This is ok, as the *order* of the GIDs in the rowmap is a better
746 // expression of the user's intent than the GIDs themselves.
747
748 Teuchos::ArrayView<const global_ordinal_type> nodeGIDs = rowMap->getLocalElementList();
749 for (size_t myRow = 0; myRow < A.getLocalNumRows(); ++myRow) {
750 local_ordinal_type local_row = myRow;
751
752 // TODO JJH 4April2014 An optimization is to use getLocalRowView. Not all matrices support this,
753 // we'd need to check via the Tpetra::RowMatrix method supportsRowViews().
754 A.getLocalRowCopy(local_row, InI, InV, NumIn); // Get Values and Indices
755
756 // Split into L and U (we don't assume that indices are ordered).
757
758 NumL = 0;
759 NumU = 0;
760 DiagFound = false;
761
762 for (size_t j = 0; j < NumIn; ++j) {
763 const local_ordinal_type k = InI[j];
764
765 if (k == local_row) {
766 DiagFound = true;
767 // Store perturbed diagonal in Tpetra::Vector D_
768 DV(local_row) += Rthresh_ * InV[j] + IFPACK2_SGN(InV[j]) * Athresh_;
769 } else if (k < 0) { // Out of range
770 TEUCHOS_TEST_FOR_EXCEPTION(
771 true, std::runtime_error,
772 "Ifpack2::RILUK::initAllValues: current "
773 "GID k = "
774 << k << " < 0. I'm not sure why this is an error; it is "
775 "probably an artifact of the undocumented assumptions of the "
776 "original implementation (likely copied and pasted from Ifpack). "
777 "Nevertheless, the code I found here insisted on this being an error "
778 "state, so I will throw an exception here.");
779 } else if (k < local_row) {
780 LI[NumL] = k;
781 LV[NumL] = InV[j];
782 NumL++;
783 } else if (Teuchos::as<size_t>(k) <= rowMap->getLocalNumElements()) {
784 UI[NumU] = k;
785 UV[NumU] = InV[j];
786 NumU++;
787 }
788 }
789
790 // Check in things for this row of L and U
791
792 if (DiagFound) {
793 ++NumNonzeroDiags;
794 } else {
795 DV(local_row) = Athresh_;
796 }
797
798 if (NumL) {
799 if (ReplaceValues) {
800 L_->replaceLocalValues(local_row, LI(0, NumL), LV(0, NumL));
801 } else {
802 // FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
803 // FIXME in this row in the column locations corresponding to UI.
804 L_->insertLocalValues(local_row, LI(0, NumL), LV(0, NumL));
805 }
806 }
807
808 if (NumU) {
809 if (ReplaceValues) {
810 U_->replaceLocalValues(local_row, UI(0, NumU), UV(0, NumU));
811 } else {
812 // FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
813 // FIXME in this row in the column locations corresponding to UI.
814 U_->insertLocalValues(local_row, UI(0, NumU), UV(0, NumU));
815 }
816 }
817 }
818
819 // At this point L and U have the values of A in the structure of L
820 // and U, and diagonal vector D
821
822 isInitialized_ = true;
823}
824
825template <class MatrixType>
826void RILUK<MatrixType>::compute_serial() {
827 // Fill L and U with numbers. This supports nonzero pattern reuse by calling
828 // initialize() once and then compute() multiple times.
829 initAllValues(*A_local_);
830
831 // MinMachNum should be officially defined, for now pick something a little
832 // bigger than IEEE underflow value
833
834 const scalar_type MinDiagonalValue = STS::rmin();
835 const scalar_type MaxDiagonalValue = STS::one() / MinDiagonalValue;
836
837 size_t NumIn, NumL, NumU;
838
839 // Get Maximum Row length
840 const size_t MaxNumEntries =
841 L_->getLocalMaxNumRowEntries() + U_->getLocalMaxNumRowEntries() + 1;
842
843 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries); // Allocate temp space
844 Teuchos::Array<scalar_type> InV(MaxNumEntries);
845 size_t num_cols = U_->getColMap()->getLocalNumElements();
846 Teuchos::Array<int> colflag(num_cols, -1);
847
848 auto DV = Kokkos::subview(D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
849
850 // Now start the factorization.
851
852 using IST = typename row_matrix_type::impl_scalar_type;
853 for (size_t i = 0; i < L_->getLocalNumRows(); ++i) {
854 local_ordinal_type local_row = i;
855 // Need some integer workspace and pointers
856 size_t NumUU;
857 local_inds_host_view_type UUI;
858 values_host_view_type UUV;
859
860 // Fill InV, InI with current row of L, D and U combined
861
862 NumIn = MaxNumEntries;
863 nonconst_local_inds_host_view_type InI_v(InI.data(), MaxNumEntries);
864 nonconst_values_host_view_type InV_v(reinterpret_cast<IST*>(InV.data()), MaxNumEntries);
865
866 L_->getLocalRowCopy(local_row, InI_v, InV_v, NumL);
867
868 InV[NumL] = DV(i); // Put in diagonal
869 InI[NumL] = local_row;
870
871 nonconst_local_inds_host_view_type InI_sub(InI.data() + NumL + 1, MaxNumEntries - NumL - 1);
872 nonconst_values_host_view_type InV_sub(reinterpret_cast<IST*>(InV.data()) + NumL + 1, MaxNumEntries - NumL - 1);
873
874 U_->getLocalRowCopy(local_row, InI_sub, InV_sub, NumU);
875 NumIn = NumL + NumU + 1;
876
877 // Set column flags
878 for (size_t j = 0; j < NumIn; ++j) {
879 colflag[InI[j]] = j;
880 }
881
882 scalar_type diagmod = STS::zero(); // Off-diagonal accumulator
883
884 for (size_t jj = 0; jj < NumL; ++jj) {
885 local_ordinal_type j = InI[jj];
886 IST multiplier = InV[jj]; // current_mults++;
887
888 InV[jj] *= static_cast<scalar_type>(DV(j));
889
890 U_->getLocalRowView(j, UUI, UUV); // View of row above
891 NumUU = UUI.size();
892
893 if (RelaxValue_ == STM::zero()) {
894 for (size_t k = 0; k < NumUU; ++k) {
895 const int kk = colflag[UUI[k]];
896 // FIXME (mfh 23 Dec 2013) Wait a second, we just set
897 // colflag above using size_t (which is generally unsigned),
898 // but now we're querying it using int (which is signed).
899 if (kk > -1) {
900 InV[kk] -= static_cast<scalar_type>(multiplier * UUV[k]);
901 }
902 }
903
904 } else {
905 for (size_t k = 0; k < NumUU; ++k) {
906 // FIXME (mfh 23 Dec 2013) Wait a second, we just set
907 // colflag above using size_t (which is generally unsigned),
908 // but now we're querying it using int (which is signed).
909 const int kk = colflag[UUI[k]];
910 if (kk > -1) {
911 InV[kk] -= static_cast<scalar_type>(multiplier * UUV[k]);
912 } else {
913 diagmod -= static_cast<scalar_type>(multiplier * UUV[k]);
914 }
915 }
916 }
917 }
918
919 if (NumL) {
920 // Replace current row of L
921 L_->replaceLocalValues(local_row, InI(0, NumL), InV(0, NumL));
922 }
923
924 DV(i) = InV[NumL]; // Extract Diagonal value
925
926 if (RelaxValue_ != STM::zero()) {
927 DV(i) += RelaxValue_ * diagmod; // Add off diagonal modifications
928 }
929
930 if (STS::magnitude(DV(i)) > STS::magnitude(MaxDiagonalValue)) {
931 if (STS::real(DV(i)) < STM::zero()) {
932 DV(i) = -MinDiagonalValue;
933 } else {
934 DV(i) = MinDiagonalValue;
935 }
936 } else {
937 DV(i) = static_cast<impl_scalar_type>(STS::one()) / DV(i); // Invert diagonal value
938 }
939
940 for (size_t j = 0; j < NumU; ++j) {
941 InV[NumL + 1 + j] *= static_cast<scalar_type>(DV(i)); // Scale U by inverse of diagonal
942 }
943
944 if (NumU) {
945 // Replace current row of L and U
946 U_->replaceLocalValues(local_row, InI(NumL + 1, NumU), InV(NumL + 1, NumU));
947 }
948
949 // Reset column flags
950 for (size_t j = 0; j < NumIn; ++j) {
951 colflag[InI[j]] = -1;
952 }
953 }
954
955 // The domain of L and the range of U are exactly their own row maps
956 // (there is no communication). The domain of U and the range of L
957 // must be the same as those of the original matrix, However if the
958 // original matrix is a VbrMatrix, these two latter maps are
959 // translation from a block map to a point map.
960 // FIXME (mfh 23 Dec 2013) Do we know that the column Map of L_ is
961 // always one-to-one?
962 L_->fillComplete(L_->getColMap(), A_local_->getRangeMap());
963 U_->fillComplete(A_local_->getDomainMap(), U_->getRowMap());
964
965 // If L_solver_ or U_solver store modified factors internally, we need to reset those
966 L_solver_->setMatrix(L_);
967 L_solver_->compute(); // NOTE: Only do compute if the pointer changed. Otherwise, do nothing
968 U_solver_->setMatrix(U_);
969 U_solver_->compute(); // NOTE: Only do compute if the pointer changed. Otherwise, do nothing
970}
971
972template <class MatrixType>
973void RILUK<MatrixType>::compute_kkspiluk() {
974 L_->resumeFill();
975 U_->resumeFill();
976
977 L_->setAllToScalar(STS::zero()); // Zero out L and U matrices
978 U_->setAllToScalar(STS::zero());
979
980 using row_map_type = typename crs_matrix_type::local_matrix_device_type::row_map_type;
981 auto lclL = L_->getLocalMatrixDevice();
982 row_map_type L_rowmap = lclL.graph.row_map;
983 auto L_entries = lclL.graph.entries;
984 auto L_values = lclL.values;
985
986 auto lclU = U_->getLocalMatrixDevice();
987 row_map_type U_rowmap = lclU.graph.row_map;
988 auto U_entries = lclU.graph.entries;
989 auto U_values = lclU.values;
990
991 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
992 KokkosSparse::spiluk_numeric(KernelHandle_.getRawPtr(), LevelOfFill_,
993 lclMtx.graph.row_map, lclMtx.graph.entries, lclMtx.values,
994 L_rowmap, L_entries, L_values, U_rowmap, U_entries, U_values);
995
996 L_->fillComplete(L_->getColMap(), A_local_->getRangeMap());
997 U_->fillComplete(A_local_->getDomainMap(), U_->getRowMap());
998
999 L_solver_->compute();
1000 U_solver_->compute();
1001}
1002
1003template <class MatrixType>
1004void RILUK<MatrixType>::compute_kkspiluk_stream() {
1005 std::vector<lno_row_view_t> L_rowmap_v(num_streams_);
1006 std::vector<lno_nonzero_view_t> L_entries_v(num_streams_);
1007 std::vector<scalar_nonzero_view_t> L_values_v(num_streams_);
1008 std::vector<lno_row_view_t> U_rowmap_v(num_streams_);
1009 std::vector<lno_nonzero_view_t> U_entries_v(num_streams_);
1010 std::vector<scalar_nonzero_view_t> U_values_v(num_streams_);
1011 std::vector<kk_handle_type*> KernelHandle_rawptr_v_(num_streams_);
1012 for (int i = 0; i < num_streams_; i++) {
1013 L_v_[i]->resumeFill();
1014 U_v_[i]->resumeFill();
1015
1016 L_v_[i]->setAllToScalar(STS::zero()); // Zero out L and U matrices
1017 U_v_[i]->setAllToScalar(STS::zero());
1018
1019 auto lclL = L_v_[i]->getLocalMatrixDevice();
1020 L_rowmap_v[i] = lclL.graph.row_map;
1021 L_entries_v[i] = lclL.graph.entries;
1022 L_values_v[i] = lclL.values;
1023
1024 auto lclU = U_v_[i]->getLocalMatrixDevice();
1025 U_rowmap_v[i] = lclU.graph.row_map;
1026 U_entries_v[i] = lclU.graph.entries;
1027 U_values_v[i] = lclU.values;
1028 KernelHandle_rawptr_v_[i] = KernelHandle_v_[i].getRawPtr();
1029 }
1030
1031 if (hasStreamReordered_) {
1032 // NOTE: For now, only do the value copying if RCM reordering is enabled in streams.
1033 // Otherwise, value copying was taken care prior to entering this function.
1034 // TODO: revisit this later.
1035 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
1036 // A_local_diagblks was already setup during initialize, just copy the corresponding
1037 // values from A_local_crs_ in parallel now.
1038 using TeamPolicy = Kokkos::TeamPolicy<execution_space>;
1039 const auto A_nrows = lclMtx.numRows();
1040 auto rows_per_block = ((A_nrows % num_streams_) == 0)
1041 ? (A_nrows / num_streams_)
1042 : (A_nrows / num_streams_ + 1);
1043 for (int i = 0; i < num_streams_; i++) {
1044 const auto start_row_offset = i * rows_per_block;
1045 auto rowptrs = A_local_diagblks_rowmap_v_[i];
1046 auto colindices = A_local_diagblks_entries_v_[i];
1047 auto values = A_local_diagblks_values_v_[i];
1048 const bool reordered = hasStreamReordered_;
1049 typename lno_nonzero_view_t::non_const_type reverse_perm = hasStreamReordered_ ? reverse_perm_v_[i] : typename lno_nonzero_view_t::non_const_type{};
1050 TeamPolicy pol(exec_space_instances_[i], A_local_diagblks_rowmap_v_[i].extent(0) - 1, Kokkos::AUTO);
1051 Kokkos::parallel_for(
1052 pol, KOKKOS_LAMBDA(const typename TeamPolicy::member_type& team) {
1053 const auto irow = team.league_rank();
1054 const auto irow_A = start_row_offset + (reordered ? reverse_perm(irow) : irow);
1055 const auto A_local_crs_row = lclMtx.rowConst(irow_A);
1056 const auto begin_row = rowptrs(irow);
1057 const auto num_entries = rowptrs(irow + 1) - begin_row;
1058 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, num_entries), [&](const int j) {
1059 const auto colidx = colindices(begin_row + j);
1060 const auto colidx_A = start_row_offset + (reordered ? reverse_perm(colidx) : colidx);
1061 // Find colidx in A_local_crs_row
1062 const int offset = KokkosSparse::findRelOffset(
1063 &A_local_crs_row.colidx(0), A_local_crs_row.length, colidx_A, 0, false);
1064 values(begin_row + j) = A_local_crs_row.value(offset);
1065 });
1066 });
1067 }
1068 }
1069
1070 KokkosSparse::Experimental::spiluk_numeric_streams(exec_space_instances_, KernelHandle_rawptr_v_, LevelOfFill_,
1071 A_local_diagblks_rowmap_v_, A_local_diagblks_entries_v_, A_local_diagblks_values_v_,
1072 L_rowmap_v, L_entries_v, L_values_v,
1073 U_rowmap_v, U_entries_v, U_values_v);
1074
1075 for (int i = 0; i < num_streams_; i++) {
1076 L_v_[i]->fillComplete();
1077 U_v_[i]->fillComplete();
1078 }
1079
1080 L_solver_->compute();
1081 U_solver_->compute();
1082}
1083
1084template <class MatrixType>
1086 using Teuchos::Array;
1087 using Teuchos::ArrayView;
1088 using Teuchos::RCP;
1089 using Teuchos::rcp;
1090 using Teuchos::rcp_const_cast;
1091 using Teuchos::rcp_dynamic_cast;
1092 const char prefix[] = "Ifpack2::RILUK::compute: ";
1093
1094 // initialize() checks this too, but it's easier for users if the
1095 // error shows them the name of the method that they actually
1096 // called, rather than the name of some internally called method.
1097 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix is null. Please "
1098 "call setMatrix() with a nonnull input before calling this method.");
1099 TEUCHOS_TEST_FOR_EXCEPTION(!A_->isFillComplete(), std::runtime_error, prefix << "The matrix is not "
1100 "fill complete. You may not invoke initialize() or compute() with this "
1101 "matrix until the matrix is fill complete. If your matrix is a "
1102 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
1103 "range Maps, if appropriate) before calling this method.");
1104
1105 if (!isInitialized()) {
1106 initialize(); // Don't count this in the compute() time
1107 }
1108
1109 Teuchos::Time timer("RILUK::compute");
1110
1111 // Start timing
1112 Teuchos::TimeMonitor timeMon(timer);
1113 double startTime = timer.wallTime();
1114
1115 isComputed_ = false;
1116
1117 if (!this->isKokkosKernelsSpiluk_) {
1118 compute_serial();
1119 } else {
1120 // Make sure values in A is picked up even in case of pattern reuse
1121 if (!A_local_crs_nc_.is_null()) {
1122 // copy entries into A_local_crs
1123 A_local_crs_nc_->resumeFill();
1124 local_ordinal_type numRows = A_local_->getLocalNumRows();
1125 nonconst_local_inds_host_view_type indices("indices", A_local_->getLocalMaxNumRowEntries());
1126 nonconst_values_host_view_type values("values", A_local_->getLocalMaxNumRowEntries());
1127 for (local_ordinal_type i = 0; i < numRows; i++) {
1128 size_t numEntries = 0;
1129 A_local_->getLocalRowCopy(i, indices, values, numEntries);
1130 A_local_crs_nc_->replaceLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()), indices.data());
1131 }
1132 A_local_crs_nc_->fillComplete(A_local_->getDomainMap(), A_local_->getRangeMap());
1133 }
1134
1135 if (!isKokkosKernelsStream_) {
1136 compute_kkspiluk();
1137 } else {
1138 // If streams are on, we potentially have to refresh A_local_diagblks_values_v_
1139 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
1140 if (!hasStreamsWithRCB_) {
1141 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, hasStreamReordered_);
1142 } else {
1143 auto A_coordinates_lcl = A_coordinates_->getLocalViewDevice(Tpetra::Access::ReadOnly);
1144 Kokkos::deep_copy(coors_rcb_, A_coordinates_lcl);
1145 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, coors_rcb_,
1146 A_local_diagblks_v_, perm_rcb_);
1147 }
1148 for (int i = 0; i < num_streams_; i++) {
1149 A_local_diagblks_values_v_[i] = A_local_diagblks_v_[i].values;
1150 }
1151
1152 compute_kkspiluk_stream();
1153 }
1154 }
1155
1156 isComputed_ = true;
1157 ++numCompute_;
1158 computeTime_ += (timer.wallTime() - startTime);
1159}
1160
1161namespace Impl {
1162template <typename MV, typename Map>
1163void resetMultiVecIfNeeded(std::unique_ptr<MV>& mv_ptr, const Map& map, const size_t numVectors, bool initialize) {
1164 if (!mv_ptr || mv_ptr->getNumVectors() != numVectors) {
1165 mv_ptr.reset(new MV(map, numVectors, initialize));
1166 }
1167}
1168} // namespace Impl
1169
1170template <class MatrixType>
1172 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
1173 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
1174 Teuchos::ETransp mode,
1175 scalar_type alpha,
1176 scalar_type beta) const {
1177 using Teuchos::RCP;
1178 using Teuchos::rcpFromRef;
1179
1180 TEUCHOS_TEST_FOR_EXCEPTION(
1181 A_.is_null(), std::runtime_error,
1182 "Ifpack2::RILUK::apply: The matrix is "
1183 "null. Please call setMatrix() with a nonnull input, then initialize() "
1184 "and compute(), before calling this method.");
1185 TEUCHOS_TEST_FOR_EXCEPTION(
1186 !isComputed(), std::runtime_error,
1187 "Ifpack2::RILUK::apply: If you have not yet called compute(), "
1188 "you must call compute() before calling this method.");
1189 TEUCHOS_TEST_FOR_EXCEPTION(
1190 X.getNumVectors() != Y.getNumVectors(), std::invalid_argument,
1191 "Ifpack2::RILUK::apply: X and Y do not have the same number of columns. "
1192 "X.getNumVectors() = "
1193 << X.getNumVectors()
1194 << " != Y.getNumVectors() = " << Y.getNumVectors() << ".");
1195 TEUCHOS_TEST_FOR_EXCEPTION(
1196 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
1197 "Ifpack2::RILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
1198 "complex Scalar type. Please talk to the Ifpack2 developers to get this "
1199 "fixed. There is a FIXME in this file about this very issue.");
1200#ifdef HAVE_IFPACK2_DEBUG
1201 {
1202 if (!isKokkosKernelsStream_) {
1203 const magnitude_type D_nrm1 = D_->norm1();
1204 TEUCHOS_TEST_FOR_EXCEPTION(STM::isnaninf(D_nrm1), std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the stored diagonal is NaN or Inf.");
1205 }
1206 Teuchos::Array<magnitude_type> norms(X.getNumVectors());
1207 X.norm1(norms());
1208 bool good = true;
1209 for (size_t j = 0; j < X.getNumVectors(); ++j) {
1210 if (STM::isnaninf(norms[j])) {
1211 good = false;
1212 break;
1213 }
1214 }
1215 TEUCHOS_TEST_FOR_EXCEPTION(!good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the input X is NaN or Inf.");
1216 }
1217#endif // HAVE_IFPACK2_DEBUG
1218
1219 const scalar_type one = STS::one();
1220 const scalar_type zero = STS::zero();
1221
1222 Teuchos::Time timer("RILUK::apply");
1223 double startTime = timer.wallTime();
1224 { // Start timing
1225 Teuchos::TimeMonitor timeMon(timer);
1226 if (alpha == one && beta == zero) {
1227 if (isKokkosKernelsSpiluk_ && isKokkosKernelsStream_ && (hasStreamReordered_ || hasStreamsWithRCB_)) {
1228 Impl::resetMultiVecIfNeeded(reordered_x_, X.getMap(), X.getNumVectors(), false);
1229 Impl::resetMultiVecIfNeeded(reordered_y_, Y.getMap(), Y.getNumVectors(), false);
1230 Kokkos::fence();
1231
1232 for (size_t j = 0; j < X.getNumVectors(); j++) {
1233 auto X_j = X.getVector(j);
1234 auto ReorderedX_j = reordered_x_->getVectorNonConst(j);
1235 auto X_lcl = X_j->getLocalViewDevice(Tpetra::Access::ReadOnly);
1236 auto ReorderedX_lcl = ReorderedX_j->getLocalViewDevice(Tpetra::Access::ReadWrite);
1237 if (hasStreamReordered_) {
1238 local_ordinal_type stream_begin = 0;
1239 local_ordinal_type stream_end;
1240 for (int i = 0; i < num_streams_; i++) {
1241 auto perm_i = perm_v_[i];
1242 stream_end = stream_begin + perm_i.extent(0);
1243 auto X_lcl_sub = Kokkos::subview(X_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1244 auto ReorderedX_lcl_sub = Kokkos::subview(ReorderedX_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1245 Kokkos::parallel_for(
1246 Kokkos::RangePolicy<execution_space>(exec_space_instances_[i], 0, static_cast<int>(perm_i.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1247 ReorderedX_lcl_sub(perm_i(ii)) = X_lcl_sub(ii);
1248 });
1249 stream_begin = stream_end;
1250 }
1251 } else { // hasStreamsWithRCB_
1252 auto perm = perm_rcb_;
1253 Kokkos::parallel_for(
1254 Kokkos::RangePolicy<execution_space>(0, static_cast<int>(X_lcl.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1255 ReorderedX_lcl(perm(ii), 0) = X_lcl(ii, 0);
1256 });
1257 }
1258 }
1259 Kokkos::fence();
1260
1261 if (mode == Teuchos::NO_TRANS) { // Solve L (U Y) = X for Y.
1262 // Solve L Y = X for Y.
1263 L_solver_->apply(*reordered_x_, Y, mode);
1264 // Solve U Y = Y for Y.
1265 U_solver_->apply(Y, *reordered_y_, mode);
1266 } else { // Solve U^P (L^P Y) = X for Y (where P is * or T).
1267 // Solve U^P Y = X for Y.
1268 U_solver_->apply(*reordered_x_, Y, mode);
1269 // Solve L^P Y = Y for Y.
1270 L_solver_->apply(Y, *reordered_y_, mode);
1271 }
1272
1273 for (size_t j = 0; j < Y.getNumVectors(); j++) {
1274 auto Y_j = Y.getVectorNonConst(j);
1275 auto ReorderedY_j = reordered_y_->getVector(j);
1276 auto Y_lcl = Y_j->getLocalViewDevice(Tpetra::Access::ReadWrite);
1277 auto ReorderedY_lcl = ReorderedY_j->getLocalViewDevice(Tpetra::Access::ReadOnly);
1278 if (hasStreamReordered_) {
1279 local_ordinal_type stream_begin = 0;
1280 local_ordinal_type stream_end;
1281 for (int i = 0; i < num_streams_; i++) {
1282 auto perm_i = perm_v_[i];
1283 stream_end = stream_begin + perm_i.extent(0);
1284 auto Y_lcl_sub = Kokkos::subview(Y_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1285 auto ReorderedY_lcl_sub = Kokkos::subview(ReorderedY_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1286 Kokkos::parallel_for(
1287 Kokkos::RangePolicy<execution_space>(exec_space_instances_[i], 0, static_cast<int>(perm_i.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1288 Y_lcl_sub(ii) = ReorderedY_lcl_sub(perm_i(ii));
1289 });
1290 stream_begin = stream_end;
1291 }
1292 } else { // hasStreamsWithRCB_
1293 auto perm = perm_rcb_;
1294 Kokkos::parallel_for(
1295 Kokkos::RangePolicy<execution_space>(0, static_cast<int>(Y_lcl.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1296 Y_lcl(ii, 0) = ReorderedY_lcl(perm(ii), 0);
1297 });
1298 }
1299 }
1300 Kokkos::fence();
1301 } else {
1302 if (mode == Teuchos::NO_TRANS) { // Solve L (D (U Y)) = X for Y.
1303#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1304 // NOTE (Nov-15-2022):
1305 // This is a workaround for Cuda >= 11.3 (using cusparseSpSV)
1306 // since cusparseSpSV_solve() does not support in-place computation
1307 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1308
1309 // Start by solving L Y_tmp = X for Y_tmp.
1310 L_solver_->apply(X, *Y_tmp_, mode);
1311
1312 if (!this->isKokkosKernelsSpiluk_) {
1313 // Solve D Y = Y. The operation lets us do this in place in Y, so we can
1314 // write "solve D Y = Y for Y."
1315 Y_tmp_->elementWiseMultiply(one, *D_, *Y_tmp_, zero);
1316 }
1317
1318 U_solver_->apply(*Y_tmp_, Y, mode); // Solve U Y = Y_tmp.
1319#else
1320 // Start by solving L Y = X for Y.
1321 L_solver_->apply(X, Y, mode);
1322
1323 if (!this->isKokkosKernelsSpiluk_) {
1324 // Solve D Y = Y. The operation lets us do this in place in Y, so we can
1325 // write "solve D Y = Y for Y."
1326 Y.elementWiseMultiply(one, *D_, Y, zero);
1327 }
1328
1329 U_solver_->apply(Y, Y, mode); // Solve U Y = Y.
1330#endif
1331 } else { // Solve U^P (D^P (L^P Y)) = X for Y (where P is * or T).
1332#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1333 // NOTE (Nov-15-2022):
1334 // This is a workaround for Cuda >= 11.3 (using cusparseSpSV)
1335 // since cusparseSpSV_solve() does not support in-place computation
1336 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1337
1338 // Start by solving U^P Y_tmp = X for Y_tmp.
1339 U_solver_->apply(X, *Y_tmp_, mode);
1340
1341 if (!this->isKokkosKernelsSpiluk_) {
1342 // Solve D^P Y = Y.
1343 //
1344 // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we
1345 // need to do an elementwise multiply with the conjugate of
1346 // D_, not just with D_ itself.
1347 Y_tmp_->elementWiseMultiply(one, *D_, *Y_tmp_, zero);
1348 }
1349
1350 L_solver_->apply(*Y_tmp_, Y, mode); // Solve L^P Y = Y_tmp.
1351#else
1352 // Start by solving U^P Y = X for Y.
1353 U_solver_->apply(X, Y, mode);
1354
1355 if (!this->isKokkosKernelsSpiluk_) {
1356 // Solve D^P Y = Y.
1357 //
1358 // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we
1359 // need to do an elementwise multiply with the conjugate of
1360 // D_, not just with D_ itself.
1361 Y.elementWiseMultiply(one, *D_, Y, zero);
1362 }
1363
1364 L_solver_->apply(Y, Y, mode); // Solve L^P Y = Y.
1365#endif
1366 }
1367 }
1368 } else { // alpha != 1 or beta != 0
1369 if (alpha == zero) {
1370 // The special case for beta == 0 ensures that if Y contains Inf
1371 // or NaN values, we replace them with 0 (following BLAS
1372 // convention), rather than multiplying them by 0 to get NaN.
1373 if (beta == zero) {
1374 Y.putScalar(zero);
1375 } else {
1376 Y.scale(beta);
1377 }
1378 } else { // alpha != zero
1379 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1380 apply(X, *Y_tmp_, mode);
1381 Y.update(alpha, *Y_tmp_, beta);
1382 }
1383 }
1384 } // end timing
1385
1386#ifdef HAVE_IFPACK2_DEBUG
1387 {
1388 Teuchos::Array<magnitude_type> norms(Y.getNumVectors());
1389 Y.norm1(norms());
1390 bool good = true;
1391 for (size_t j = 0; j < Y.getNumVectors(); ++j) {
1392 if (STM::isnaninf(norms[j])) {
1393 good = false;
1394 break;
1395 }
1396 }
1397 TEUCHOS_TEST_FOR_EXCEPTION(!good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the output Y is NaN or Inf.");
1398 }
1399#endif // HAVE_IFPACK2_DEBUG
1400
1401 ++numApply_;
1402 applyTime_ += (timer.wallTime() - startTime);
1403}
1404
1405// VINH comment out since multiply() is not needed anywhere
1406// template<class MatrixType>
1407// void RILUK<MatrixType>::
1408// multiply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1409// Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1410// const Teuchos::ETransp mode) const
1411//{
1412// const scalar_type zero = STS::zero ();
1413// const scalar_type one = STS::one ();
1414//
1415// if (mode != Teuchos::NO_TRANS) {
1416// U_->apply (X, Y, mode); //
1417// Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1418//
1419// // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we need
1420// // to do an elementwise multiply with the conjugate of D_, not
1421// // just with D_ itself.
1422// Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1423//
1424// MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y
1425// L_->apply (Y_tmp, Y, mode);
1426// Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1427// }
1428// else {
1429// L_->apply (X, Y, mode);
1430// Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1431// Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1432// MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y1
1433// U_->apply (Y_tmp, Y, mode);
1434// Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1435// }
1436// }
1437
1438template <class MatrixType>
1440 std::ostringstream os;
1441
1442 // Output is a valid YAML dictionary in flow style. If you don't
1443 // like everything on a single line, you should call describe()
1444 // instead.
1445 os << "\"Ifpack2::RILUK\": {";
1446 os << "Initialized: " << (isInitialized() ? "true" : "false") << ", "
1447 << "Computed: " << (isComputed() ? "true" : "false") << ", ";
1448
1449 os << "Level-of-fill: " << getLevelOfFill() << ", ";
1450
1451 if (isKokkosKernelsSpiluk_) os << "KK-SPILUK, ";
1452 if (isKokkosKernelsStream_) os << "KK-Stream, ";
1453
1454 if (A_.is_null()) {
1455 os << "Matrix: null";
1456 } else {
1457 os << "Global matrix dimensions: ["
1458 << A_->getGlobalNumRows() << ", " << A_->getGlobalNumCols() << "]"
1459 << ", Global nnz: " << A_->getGlobalNumEntries();
1460 }
1461
1462 if (!L_solver_.is_null()) os << ", " << L_solver_->description();
1463 if (!U_solver_.is_null()) os << ", " << U_solver_->description();
1464
1465 os << "}";
1466 return os.str();
1467}
1468
1469} // namespace Ifpack2
1470
1471#define IFPACK2_RILUK_INSTANT(S, LO, GO, N) \
1472 template class Ifpack2::RILUK<Tpetra::RowMatrix<S, LO, GO, N> >;
1473
1474#endif
Declaration of MDF interface.
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition Ifpack2_IlukGraph.hpp:67
Access only local rows and columns of a sparse matrix.
Definition Ifpack2_LocalFilter_decl.hpp:128
"Preconditioner" that solves local sparse triangular systems.
Definition Ifpack2_LocalSparseTriangularSolver_decl.hpp:54
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition Ifpack2_RILUK_decl.hpp:218
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:468
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition Ifpack2_RILUK_decl.hpp:255
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition Ifpack2_RILUK_def.hpp:99
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
Definition Ifpack2_RILUK_def.hpp:1172
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_RILUK_decl.hpp:233
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_RILUK_def.hpp:240
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:191
std::string description() const
A one-line description of this object.
Definition Ifpack2_RILUK_def.hpp:1439
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_RILUK_def.hpp:122
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:221
Teuchos::RCP< const coord_type > getCoord() const
Get the coordinates associated with the input matrix's rows.
Definition Ifpack2_RILUK_def.hpp:432
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > & getD() const
Return the diagonal entries of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:178
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:230
void setCoord(const Teuchos::RCP< const coord_type > &A_coordinates)
Set the matrix rows' coordinates.
Definition Ifpack2_RILUK_def.hpp:154
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_RILUK_def.hpp:420
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition Ifpack2_RILUK_def.hpp:426
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_RILUK_def.hpp:220
void setParameters(const Teuchos::ParameterList &params)
Definition Ifpack2_RILUK_def.hpp:299
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:224
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:227
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:162
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:239
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_RILUK_def.hpp:203
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:1085
static Teuchos::RCP< const row_matrix_type > makeLocalFilter(const Teuchos::RCP< const row_matrix_type > &A)
Return A, wrapped in a LocalFilter, if necessary.
Definition Ifpack2_RILUK_def.hpp:438
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40