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#if KOKKOS_VERSION >= 50100
574 reverse_perm_rcb_ = perm_view_t(Kokkos::view_alloc(Kokkos::WithoutInitializing, "reverse_perm_rcb_"), A_coordinates_lcl.extent(0));
575#endif
576 coors_rcb_ = coors_view_t(Kokkos::view_alloc(Kokkos::WithoutInitializing, "coors_rcb_"), A_coordinates_lcl.extent(0), A_coordinates_lcl.extent(1));
577 Kokkos::deep_copy(coors_rcb_, A_coordinates_lcl);
578#if KOKKOS_VERSION >= 50100
579 local_ordinal_type n_levels = static_cast<local_ordinal_type>(std::log2(static_cast<double>(num_streams_)) + 1);
580 partition_sizes_rcb_ = KokkosGraph::Experimental::recursive_coordinate_bisection(coors_rcb_, perm_rcb_, reverse_perm_rcb_, n_levels);
581 KokkosSparse::Experimental::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, perm_rcb_, reverse_perm_rcb_, partition_sizes_rcb_,
582 A_local_diagblks_v_);
583#else
584 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, coors_rcb_,
585 A_local_diagblks_v_, perm_rcb_);
586#endif
587 } else {
588#if KOKKOS_VERSION >= 50100
589 KokkosSparse::Experimental::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_);
590#else
591 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_);
592#endif
593 }
594 } else {
595#if KOKKOS_VERSION >= 50100
596 perm_v_ = KokkosSparse::Experimental::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, true);
597#else
598 perm_v_ = KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, true);
599#endif
600 reverse_perm_v_.resize(perm_v_.size());
601 for (size_t istream = 0; istream < perm_v_.size(); ++istream) {
602 using perm_type = typename lno_nonzero_view_t::non_const_type;
603 const auto perm = perm_v_[istream];
604 const auto perm_length = perm.extent(0);
605 perm_type reverse_perm(
606 Kokkos::view_alloc(Kokkos::WithoutInitializing, "reverse_perm"),
607 perm_length);
608 Kokkos::parallel_for(
609 Kokkos::RangePolicy<execution_space>(exec_space_instances_[istream], 0, perm_length),
610 KOKKOS_LAMBDA(const local_ordinal_type ii) {
611 reverse_perm(perm(ii)) = ii;
612 });
613 reverse_perm_v_[istream] = reverse_perm;
614 }
615 }
616
617 A_local_diagblks_rowmap_v_ = std::vector<lno_row_view_t>(num_streams_);
618 A_local_diagblks_entries_v_ = std::vector<lno_nonzero_view_t>(num_streams_);
619 A_local_diagblks_values_v_ = std::vector<scalar_nonzero_view_t>(num_streams_);
620
621 for (int i = 0; i < num_streams_; i++) {
622 A_local_diagblks_rowmap_v_[i] = A_local_diagblks_v_[i].graph.row_map;
623 A_local_diagblks_entries_v_[i] = A_local_diagblks_v_[i].graph.entries;
624 A_local_diagblks_values_v_[i] = A_local_diagblks_v_[i].values;
625
626 Teuchos::RCP<const crs_map_type> A_local_diagblks_RowMap = rcp(new crs_map_type(A_local_diagblks_v_[i].numRows(),
627 A_local_diagblks_v_[i].numRows(),
628 A_local_crs_->getRowMap()->getComm()));
629 Teuchos::RCP<const crs_map_type> A_local_diagblks_ColMap = rcp(new crs_map_type(A_local_diagblks_v_[i].numCols(),
630 A_local_diagblks_v_[i].numCols(),
631 A_local_crs_->getColMap()->getComm()));
632 Teuchos::RCP<crs_matrix_type> A_local_diagblks = rcp(new crs_matrix_type(A_local_diagblks_RowMap,
633 A_local_diagblks_ColMap,
634 A_local_diagblks_v_[i]));
635 Graph_v_[i] = rcp(new Ifpack2::IlukGraph<crs_graph_type, kk_handle_type>(A_local_diagblks->getCrsGraph(),
636 LevelOfFill_, 0, Overalloc_));
637 }
638 }
639 }
640
641 if (this->isKokkosKernelsSpiluk_) {
642 if (!isKokkosKernelsStream_) {
643 this->KernelHandle_ = Teuchos::rcp(new kk_handle_type());
644 KernelHandle_->create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1,
645 A_local_->getLocalNumRows(),
646 2 * A_local_->getLocalNumEntries() * (LevelOfFill_ + 1),
647 2 * A_local_->getLocalNumEntries() * (LevelOfFill_ + 1));
648 Graph_->initialize(KernelHandle_); // this calls spiluk_symbolic
649 } else {
650 KernelHandle_v_ = std::vector<Teuchos::RCP<kk_handle_type> >(num_streams_);
651 for (int i = 0; i < num_streams_; i++) {
652 KernelHandle_v_[i] = Teuchos::rcp(new kk_handle_type());
653 KernelHandle_v_[i]->create_spiluk_handle(KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1,
654 A_local_diagblks_v_[i].numRows(),
655 2 * A_local_diagblks_v_[i].nnz() * (LevelOfFill_ + 1),
656 2 * A_local_diagblks_v_[i].nnz() * (LevelOfFill_ + 1));
657 Graph_v_[i]->initialize(KernelHandle_v_[i]); // this calls spiluk_symbolic
658 }
659 }
660 } else {
661 Graph_->initialize();
662 }
663
664 allocate_L_and_U();
665 checkOrderingConsistency(*A_local_);
666 if (!isKokkosKernelsStream_) {
667 L_solver_->setMatrix(L_);
668 } else {
669 L_solver_->setStreamInfo(isKokkosKernelsStream_, num_streams_, exec_space_instances_);
670 L_solver_->setMatrices(L_v_);
671 }
672 L_solver_->initialize();
673
674 if (!isKokkosKernelsStream_) {
675 U_solver_->setMatrix(U_);
676 } else {
677 U_solver_->setStreamInfo(isKokkosKernelsStream_, num_streams_, exec_space_instances_);
678 U_solver_->setMatrices(U_v_);
679 }
680 U_solver_->initialize();
681
682 // Do not call initAllValues. compute() always calls initAllValues to
683 // fill L and U with possibly new numbers. initialize() is concerned
684 // only with the nonzero pattern.
685 } // Stop timing
686
687 isInitialized_ = true;
688 ++numInitialize_;
689 initializeTime_ += (timer.wallTime() - startTime);
690}
691
692template <class MatrixType>
693void RILUK<MatrixType>::
694 checkOrderingConsistency(const row_matrix_type& A) {
695 // First check that the local row map ordering is the same as the local portion of the column map.
696 // The extraction of the strictly lower/upper parts of A, as well as the factorization,
697 // implicitly assume that this is the case.
698 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getLocalElementList();
699 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getLocalElementList();
700 bool gidsAreConsistentlyOrdered = true;
701 global_ordinal_type indexOfInconsistentGID = 0;
702 for (global_ordinal_type i = 0; i < rowGIDs.size(); ++i) {
703 if (rowGIDs[i] != colGIDs[i]) {
704 gidsAreConsistentlyOrdered = false;
705 indexOfInconsistentGID = i;
706 break;
707 }
708 }
709 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered == false, std::runtime_error,
710 "The ordering of the local GIDs in the row and column maps is not the same"
711 << std::endl
712 << "at index " << indexOfInconsistentGID
713 << ". Consistency is required, as all calculations are done with"
714 << std::endl
715 << "local indexing.");
716}
717
718template <class MatrixType>
719void RILUK<MatrixType>::
720 initAllValues(const row_matrix_type& A) {
721 using Teuchos::ArrayRCP;
722 using Teuchos::Comm;
723 using Teuchos::ptr;
724 using Teuchos::RCP;
725 using Teuchos::REDUCE_SUM;
726 using Teuchos::reduceAll;
727 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
728
729 size_t NumIn = 0, NumL = 0, NumU = 0;
730 bool DiagFound = false;
731 size_t NumNonzeroDiags = 0;
732 size_t MaxNumEntries = A.getGlobalMaxNumRowEntries();
733
734 // Allocate temporary space for extracting the strictly
735 // lower and upper parts of the matrix A.
736 nonconst_local_inds_host_view_type InI("InI", MaxNumEntries);
737 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
738 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
739 nonconst_values_host_view_type InV("InV", MaxNumEntries);
740 Teuchos::Array<scalar_type> LV(MaxNumEntries);
741 Teuchos::Array<scalar_type> UV(MaxNumEntries);
742
743 // Check if values should be inserted or replaced
744 const bool ReplaceValues = L_->isStaticGraph() || L_->isLocallyIndexed();
745
746 L_->resumeFill();
747 U_->resumeFill();
748 if (ReplaceValues) {
749 L_->setAllToScalar(STS::zero()); // Zero out L and U matrices
750 U_->setAllToScalar(STS::zero());
751 }
752
753 D_->putScalar(STS::zero()); // Set diagonal values to zero
754 auto DV = Kokkos::subview(D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
755
756 RCP<const map_type> rowMap = L_->getRowMap();
757
758 // First we copy the user's matrix into L and U, regardless of fill level.
759 // It is important to note that L and U are populated using local indices.
760 // This means that if the row map GIDs are not monotonically increasing
761 // (i.e., permuted or gappy), then the strictly lower (upper) part of the
762 // matrix is not the one that you would get if you based L (U) on GIDs.
763 // This is ok, as the *order* of the GIDs in the rowmap is a better
764 // expression of the user's intent than the GIDs themselves.
765
766 Teuchos::ArrayView<const global_ordinal_type> nodeGIDs = rowMap->getLocalElementList();
767 for (size_t myRow = 0; myRow < A.getLocalNumRows(); ++myRow) {
768 local_ordinal_type local_row = myRow;
769
770 // TODO JJH 4April2014 An optimization is to use getLocalRowView. Not all matrices support this,
771 // we'd need to check via the Tpetra::RowMatrix method supportsRowViews().
772 A.getLocalRowCopy(local_row, InI, InV, NumIn); // Get Values and Indices
773
774 // Split into L and U (we don't assume that indices are ordered).
775
776 NumL = 0;
777 NumU = 0;
778 DiagFound = false;
779
780 for (size_t j = 0; j < NumIn; ++j) {
781 const local_ordinal_type k = InI[j];
782
783 if (k == local_row) {
784 DiagFound = true;
785 // Store perturbed diagonal in Tpetra::Vector D_
786 DV(local_row) += Rthresh_ * InV[j] + IFPACK2_SGN(InV[j]) * Athresh_;
787 } else if (k < 0) { // Out of range
788 TEUCHOS_TEST_FOR_EXCEPTION(
789 true, std::runtime_error,
790 "Ifpack2::RILUK::initAllValues: current "
791 "GID k = "
792 << k << " < 0. I'm not sure why this is an error; it is "
793 "probably an artifact of the undocumented assumptions of the "
794 "original implementation (likely copied and pasted from Ifpack). "
795 "Nevertheless, the code I found here insisted on this being an error "
796 "state, so I will throw an exception here.");
797 } else if (k < local_row) {
798 LI[NumL] = k;
799 LV[NumL] = InV[j];
800 NumL++;
801 } else if (Teuchos::as<size_t>(k) <= rowMap->getLocalNumElements()) {
802 UI[NumU] = k;
803 UV[NumU] = InV[j];
804 NumU++;
805 }
806 }
807
808 // Check in things for this row of L and U
809
810 if (DiagFound) {
811 ++NumNonzeroDiags;
812 } else {
813 DV(local_row) = Athresh_;
814 }
815
816 if (NumL) {
817 if (ReplaceValues) {
818 L_->replaceLocalValues(local_row, LI(0, NumL), LV(0, NumL));
819 } else {
820 // FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
821 // FIXME in this row in the column locations corresponding to UI.
822 L_->insertLocalValues(local_row, LI(0, NumL), LV(0, NumL));
823 }
824 }
825
826 if (NumU) {
827 if (ReplaceValues) {
828 U_->replaceLocalValues(local_row, UI(0, NumU), UV(0, NumU));
829 } else {
830 // FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
831 // FIXME in this row in the column locations corresponding to UI.
832 U_->insertLocalValues(local_row, UI(0, NumU), UV(0, NumU));
833 }
834 }
835 }
836
837 // At this point L and U have the values of A in the structure of L
838 // and U, and diagonal vector D
839
840 isInitialized_ = true;
841}
842
843template <class MatrixType>
844void RILUK<MatrixType>::compute_serial() {
845 // Fill L and U with numbers. This supports nonzero pattern reuse by calling
846 // initialize() once and then compute() multiple times.
847 initAllValues(*A_local_);
848
849 // MinMachNum should be officially defined, for now pick something a little
850 // bigger than IEEE underflow value
851
852 const scalar_type MinDiagonalValue = STS::rmin();
853 const scalar_type MaxDiagonalValue = STS::one() / MinDiagonalValue;
854
855 size_t NumIn, NumL, NumU;
856
857 // Get Maximum Row length
858 const size_t MaxNumEntries =
859 L_->getLocalMaxNumRowEntries() + U_->getLocalMaxNumRowEntries() + 1;
860
861 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries); // Allocate temp space
862 Teuchos::Array<scalar_type> InV(MaxNumEntries);
863 size_t num_cols = U_->getColMap()->getLocalNumElements();
864 Teuchos::Array<int> colflag(num_cols, -1);
865
866 auto DV = Kokkos::subview(D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
867
868 // Now start the factorization.
869
870 using IST = typename row_matrix_type::impl_scalar_type;
871 for (size_t i = 0; i < L_->getLocalNumRows(); ++i) {
872 local_ordinal_type local_row = i;
873 // Need some integer workspace and pointers
874 size_t NumUU;
875 local_inds_host_view_type UUI;
876 values_host_view_type UUV;
877
878 // Fill InV, InI with current row of L, D and U combined
879
880 NumIn = MaxNumEntries;
881 nonconst_local_inds_host_view_type InI_v(InI.data(), MaxNumEntries);
882 nonconst_values_host_view_type InV_v(reinterpret_cast<IST*>(InV.data()), MaxNumEntries);
883
884 L_->getLocalRowCopy(local_row, InI_v, InV_v, NumL);
885
886 InV[NumL] = DV(i); // Put in diagonal
887 InI[NumL] = local_row;
888
889 nonconst_local_inds_host_view_type InI_sub(InI.data() + NumL + 1, MaxNumEntries - NumL - 1);
890 nonconst_values_host_view_type InV_sub(reinterpret_cast<IST*>(InV.data()) + NumL + 1, MaxNumEntries - NumL - 1);
891
892 U_->getLocalRowCopy(local_row, InI_sub, InV_sub, NumU);
893 NumIn = NumL + NumU + 1;
894
895 // Set column flags
896 for (size_t j = 0; j < NumIn; ++j) {
897 colflag[InI[j]] = j;
898 }
899
900 scalar_type diagmod = STS::zero(); // Off-diagonal accumulator
901
902 for (size_t jj = 0; jj < NumL; ++jj) {
903 local_ordinal_type j = InI[jj];
904 IST multiplier = InV[jj]; // current_mults++;
905
906 InV[jj] *= static_cast<scalar_type>(DV(j));
907
908 U_->getLocalRowView(j, UUI, UUV); // View of row above
909 NumUU = UUI.size();
910
911 if (RelaxValue_ == STM::zero()) {
912 for (size_t k = 0; k < NumUU; ++k) {
913 const int kk = colflag[UUI[k]];
914 // FIXME (mfh 23 Dec 2013) Wait a second, we just set
915 // colflag above using size_t (which is generally unsigned),
916 // but now we're querying it using int (which is signed).
917 if (kk > -1) {
918 InV[kk] -= static_cast<scalar_type>(multiplier * UUV[k]);
919 }
920 }
921
922 } else {
923 for (size_t k = 0; k < NumUU; ++k) {
924 // FIXME (mfh 23 Dec 2013) Wait a second, we just set
925 // colflag above using size_t (which is generally unsigned),
926 // but now we're querying it using int (which is signed).
927 const int kk = colflag[UUI[k]];
928 if (kk > -1) {
929 InV[kk] -= static_cast<scalar_type>(multiplier * UUV[k]);
930 } else {
931 diagmod -= static_cast<scalar_type>(multiplier * UUV[k]);
932 }
933 }
934 }
935 }
936
937 if (NumL) {
938 // Replace current row of L
939 L_->replaceLocalValues(local_row, InI(0, NumL), InV(0, NumL));
940 }
941
942 DV(i) = InV[NumL]; // Extract Diagonal value
943
944 if (RelaxValue_ != STM::zero()) {
945 DV(i) += RelaxValue_ * diagmod; // Add off diagonal modifications
946 }
947
948 if (STS::magnitude(DV(i)) > STS::magnitude(MaxDiagonalValue)) {
949 if (STS::real(DV(i)) < STM::zero()) {
950 DV(i) = -MinDiagonalValue;
951 } else {
952 DV(i) = MinDiagonalValue;
953 }
954 } else {
955 DV(i) = static_cast<impl_scalar_type>(STS::one()) / DV(i); // Invert diagonal value
956 }
957
958 for (size_t j = 0; j < NumU; ++j) {
959 InV[NumL + 1 + j] *= static_cast<scalar_type>(DV(i)); // Scale U by inverse of diagonal
960 }
961
962 if (NumU) {
963 // Replace current row of L and U
964 U_->replaceLocalValues(local_row, InI(NumL + 1, NumU), InV(NumL + 1, NumU));
965 }
966
967 // Reset column flags
968 for (size_t j = 0; j < NumIn; ++j) {
969 colflag[InI[j]] = -1;
970 }
971 }
972
973 // The domain of L and the range of U are exactly their own row maps
974 // (there is no communication). The domain of U and the range of L
975 // must be the same as those of the original matrix, However if the
976 // original matrix is a VbrMatrix, these two latter maps are
977 // translation from a block map to a point map.
978 // FIXME (mfh 23 Dec 2013) Do we know that the column Map of L_ is
979 // always one-to-one?
980 L_->fillComplete(L_->getColMap(), A_local_->getRangeMap());
981 U_->fillComplete(A_local_->getDomainMap(), U_->getRowMap());
982
983 // If L_solver_ or U_solver store modified factors internally, we need to reset those
984 L_solver_->setMatrix(L_);
985 L_solver_->compute(); // NOTE: Only do compute if the pointer changed. Otherwise, do nothing
986 U_solver_->setMatrix(U_);
987 U_solver_->compute(); // NOTE: Only do compute if the pointer changed. Otherwise, do nothing
988}
989
990template <class MatrixType>
991void RILUK<MatrixType>::compute_kkspiluk() {
992 L_->resumeFill();
993 U_->resumeFill();
994
995 L_->setAllToScalar(STS::zero()); // Zero out L and U matrices
996 U_->setAllToScalar(STS::zero());
997
998 using row_map_type = typename crs_matrix_type::local_matrix_device_type::row_map_type;
999 auto lclL = L_->getLocalMatrixDevice();
1000 row_map_type L_rowmap = lclL.graph.row_map;
1001 auto L_entries = lclL.graph.entries;
1002 auto L_values = lclL.values;
1003
1004 auto lclU = U_->getLocalMatrixDevice();
1005 row_map_type U_rowmap = lclU.graph.row_map;
1006 auto U_entries = lclU.graph.entries;
1007 auto U_values = lclU.values;
1008
1009 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
1010 KokkosSparse::spiluk_numeric(KernelHandle_.getRawPtr(), LevelOfFill_,
1011 lclMtx.graph.row_map, lclMtx.graph.entries, lclMtx.values,
1012 L_rowmap, L_entries, L_values, U_rowmap, U_entries, U_values);
1013
1014 L_->fillComplete(L_->getColMap(), A_local_->getRangeMap());
1015 U_->fillComplete(A_local_->getDomainMap(), U_->getRowMap());
1016
1017 L_solver_->compute();
1018 U_solver_->compute();
1019}
1020
1021template <class MatrixType>
1022void RILUK<MatrixType>::compute_kkspiluk_stream() {
1023 std::vector<lno_row_view_t> L_rowmap_v(num_streams_);
1024 std::vector<lno_nonzero_view_t> L_entries_v(num_streams_);
1025 std::vector<scalar_nonzero_view_t> L_values_v(num_streams_);
1026 std::vector<lno_row_view_t> U_rowmap_v(num_streams_);
1027 std::vector<lno_nonzero_view_t> U_entries_v(num_streams_);
1028 std::vector<scalar_nonzero_view_t> U_values_v(num_streams_);
1029 std::vector<kk_handle_type*> KernelHandle_rawptr_v_(num_streams_);
1030 for (int i = 0; i < num_streams_; i++) {
1031 L_v_[i]->resumeFill();
1032 U_v_[i]->resumeFill();
1033
1034 L_v_[i]->setAllToScalar(STS::zero()); // Zero out L and U matrices
1035 U_v_[i]->setAllToScalar(STS::zero());
1036
1037 auto lclL = L_v_[i]->getLocalMatrixDevice();
1038 L_rowmap_v[i] = lclL.graph.row_map;
1039 L_entries_v[i] = lclL.graph.entries;
1040 L_values_v[i] = lclL.values;
1041
1042 auto lclU = U_v_[i]->getLocalMatrixDevice();
1043 U_rowmap_v[i] = lclU.graph.row_map;
1044 U_entries_v[i] = lclU.graph.entries;
1045 U_values_v[i] = lclU.values;
1046 KernelHandle_rawptr_v_[i] = KernelHandle_v_[i].getRawPtr();
1047 }
1048
1049 if (hasStreamReordered_) {
1050 // NOTE: For now, only do the value copying if RCM reordering is enabled in streams.
1051 // Otherwise, value copying was taken care prior to entering this function.
1052 // TODO: revisit this later.
1053 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
1054 // A_local_diagblks was already setup during initialize, just copy the corresponding
1055 // values from A_local_crs_ in parallel now.
1056 using TeamPolicy = Kokkos::TeamPolicy<execution_space>;
1057 const auto A_nrows = lclMtx.numRows();
1058 auto rows_per_block = ((A_nrows % num_streams_) == 0)
1059 ? (A_nrows / num_streams_)
1060 : (A_nrows / num_streams_ + 1);
1061 for (int i = 0; i < num_streams_; i++) {
1062 const auto start_row_offset = i * rows_per_block;
1063 auto rowptrs = A_local_diagblks_rowmap_v_[i];
1064 auto colindices = A_local_diagblks_entries_v_[i];
1065 auto values = A_local_diagblks_values_v_[i];
1066 const bool reordered = hasStreamReordered_;
1067 typename lno_nonzero_view_t::non_const_type reverse_perm = hasStreamReordered_ ? reverse_perm_v_[i] : typename lno_nonzero_view_t::non_const_type{};
1068 TeamPolicy pol(exec_space_instances_[i], A_local_diagblks_rowmap_v_[i].extent(0) - 1, Kokkos::AUTO);
1069 Kokkos::parallel_for(
1070 pol, KOKKOS_LAMBDA(const typename TeamPolicy::member_type& team) {
1071 const auto irow = team.league_rank();
1072 const auto irow_A = start_row_offset + (reordered ? reverse_perm(irow) : irow);
1073 const auto A_local_crs_row = lclMtx.rowConst(irow_A);
1074 const auto begin_row = rowptrs(irow);
1075 const auto num_entries = rowptrs(irow + 1) - begin_row;
1076 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, num_entries), [&](const int j) {
1077 const auto colidx = colindices(begin_row + j);
1078 const auto colidx_A = start_row_offset + (reordered ? reverse_perm(colidx) : colidx);
1079 // Find colidx in A_local_crs_row
1080 const int offset = KokkosSparse::findRelOffset(
1081 &A_local_crs_row.colidx(0), A_local_crs_row.length, colidx_A, 0, false);
1082 values(begin_row + j) = A_local_crs_row.value(offset);
1083 });
1084 });
1085 }
1086 }
1087
1088 KokkosSparse::Experimental::spiluk_numeric_streams(exec_space_instances_, KernelHandle_rawptr_v_, LevelOfFill_,
1089 A_local_diagblks_rowmap_v_, A_local_diagblks_entries_v_, A_local_diagblks_values_v_,
1090 L_rowmap_v, L_entries_v, L_values_v,
1091 U_rowmap_v, U_entries_v, U_values_v);
1092
1093 for (int i = 0; i < num_streams_; i++) {
1094 L_v_[i]->fillComplete();
1095 U_v_[i]->fillComplete();
1096 }
1097
1098 L_solver_->compute();
1099 U_solver_->compute();
1100}
1101
1102template <class MatrixType>
1104 using Teuchos::Array;
1105 using Teuchos::ArrayView;
1106 using Teuchos::RCP;
1107 using Teuchos::rcp;
1108 using Teuchos::rcp_const_cast;
1109 using Teuchos::rcp_dynamic_cast;
1110 const char prefix[] = "Ifpack2::RILUK::compute: ";
1111
1112 // initialize() checks this too, but it's easier for users if the
1113 // error shows them the name of the method that they actually
1114 // called, rather than the name of some internally called method.
1115 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix is null. Please "
1116 "call setMatrix() with a nonnull input before calling this method.");
1117 TEUCHOS_TEST_FOR_EXCEPTION(!A_->isFillComplete(), std::runtime_error, prefix << "The matrix is not "
1118 "fill complete. You may not invoke initialize() or compute() with this "
1119 "matrix until the matrix is fill complete. If your matrix is a "
1120 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
1121 "range Maps, if appropriate) before calling this method.");
1122
1123 if (!isInitialized()) {
1124 initialize(); // Don't count this in the compute() time
1125 }
1126
1127 Teuchos::Time timer("RILUK::compute");
1128
1129 // Start timing
1130 Teuchos::TimeMonitor timeMon(timer);
1131 double startTime = timer.wallTime();
1132
1133 isComputed_ = false;
1134
1135 if (!this->isKokkosKernelsSpiluk_) {
1136 compute_serial();
1137 } else {
1138 // Make sure values in A is picked up even in case of pattern reuse
1139 if (!A_local_crs_nc_.is_null()) {
1140 // copy entries into A_local_crs
1141 A_local_crs_nc_->resumeFill();
1142 local_ordinal_type numRows = A_local_->getLocalNumRows();
1143 nonconst_local_inds_host_view_type indices("indices", A_local_->getLocalMaxNumRowEntries());
1144 nonconst_values_host_view_type values("values", A_local_->getLocalMaxNumRowEntries());
1145 for (local_ordinal_type i = 0; i < numRows; i++) {
1146 size_t numEntries = 0;
1147 A_local_->getLocalRowCopy(i, indices, values, numEntries);
1148 A_local_crs_nc_->replaceLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()), indices.data());
1149 }
1150 A_local_crs_nc_->fillComplete(A_local_->getDomainMap(), A_local_->getRangeMap());
1151 }
1152
1153 if (!isKokkosKernelsStream_) {
1154 compute_kkspiluk();
1155 } else {
1156 // If streams are on, we potentially have to refresh A_local_diagblks_values_v_
1157 auto lclMtx = A_local_crs_->getLocalMatrixDevice();
1158 if (!hasStreamsWithRCB_) {
1159#if KOKKOS_VERSION >= 50100
1160 KokkosSparse::Experimental::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, hasStreamReordered_);
1161#else
1162 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_sequential(lclMtx, A_local_diagblks_v_, hasStreamReordered_);
1163#endif
1164 } else {
1165#if KOKKOS_VERSION >= 50100
1166 KokkosSparse::Experimental::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, perm_rcb_, reverse_perm_rcb_, partition_sizes_rcb_,
1167 A_local_diagblks_v_);
1168#else
1169 auto A_coordinates_lcl = A_coordinates_->getLocalViewDevice(Tpetra::Access::ReadOnly);
1170 Kokkos::deep_copy(coors_rcb_, A_coordinates_lcl);
1171 KokkosSparse::Impl::kk_extract_diagonal_blocks_crsmatrix_with_rcb_sequential(lclMtx, coors_rcb_,
1172 A_local_diagblks_v_, perm_rcb_);
1173#endif
1174 }
1175 for (int i = 0; i < num_streams_; i++) {
1176 A_local_diagblks_values_v_[i] = A_local_diagblks_v_[i].values;
1177 }
1178
1179 compute_kkspiluk_stream();
1180 }
1181 }
1182
1183 isComputed_ = true;
1184 ++numCompute_;
1185 computeTime_ += (timer.wallTime() - startTime);
1186}
1187
1188namespace Impl {
1189template <typename MV, typename Map>
1190void resetMultiVecIfNeeded(std::unique_ptr<MV>& mv_ptr, const Map& map, const size_t numVectors, bool initialize) {
1191 if (!mv_ptr || mv_ptr->getNumVectors() != numVectors) {
1192 mv_ptr.reset(new MV(map, numVectors, initialize));
1193 }
1194}
1195} // namespace Impl
1196
1197template <class MatrixType>
1199 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
1200 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
1201 Teuchos::ETransp mode,
1202 scalar_type alpha,
1203 scalar_type beta) const {
1204 using Teuchos::RCP;
1205 using Teuchos::rcpFromRef;
1206
1207 TEUCHOS_TEST_FOR_EXCEPTION(
1208 A_.is_null(), std::runtime_error,
1209 "Ifpack2::RILUK::apply: The matrix is "
1210 "null. Please call setMatrix() with a nonnull input, then initialize() "
1211 "and compute(), before calling this method.");
1212 TEUCHOS_TEST_FOR_EXCEPTION(
1213 !isComputed(), std::runtime_error,
1214 "Ifpack2::RILUK::apply: If you have not yet called compute(), "
1215 "you must call compute() before calling this method.");
1216 TEUCHOS_TEST_FOR_EXCEPTION(
1217 X.getNumVectors() != Y.getNumVectors(), std::invalid_argument,
1218 "Ifpack2::RILUK::apply: X and Y do not have the same number of columns. "
1219 "X.getNumVectors() = "
1220 << X.getNumVectors()
1221 << " != Y.getNumVectors() = " << Y.getNumVectors() << ".");
1222 TEUCHOS_TEST_FOR_EXCEPTION(
1223 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
1224 "Ifpack2::RILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
1225 "complex Scalar type. Please talk to the Ifpack2 developers to get this "
1226 "fixed. There is a FIXME in this file about this very issue.");
1227#ifdef HAVE_IFPACK2_DEBUG
1228 {
1229 if (!isKokkosKernelsStream_) {
1230 const magnitude_type D_nrm1 = D_->norm1();
1231 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.");
1232 }
1233 Teuchos::Array<magnitude_type> norms(X.getNumVectors());
1234 X.norm1(norms());
1235 bool good = true;
1236 for (size_t j = 0; j < X.getNumVectors(); ++j) {
1237 if (STM::isnaninf(norms[j])) {
1238 good = false;
1239 break;
1240 }
1241 }
1242 TEUCHOS_TEST_FOR_EXCEPTION(!good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the input X is NaN or Inf.");
1243 }
1244#endif // HAVE_IFPACK2_DEBUG
1245
1246 const scalar_type one = STS::one();
1247 const scalar_type zero = STS::zero();
1248
1249 Teuchos::Time timer("RILUK::apply");
1250 double startTime = timer.wallTime();
1251 { // Start timing
1252 Teuchos::TimeMonitor timeMon(timer);
1253 if (alpha == one && beta == zero) {
1254 if (isKokkosKernelsSpiluk_ && isKokkosKernelsStream_ && (hasStreamReordered_ || hasStreamsWithRCB_)) {
1255 Impl::resetMultiVecIfNeeded(reordered_x_, X.getMap(), X.getNumVectors(), false);
1256 Impl::resetMultiVecIfNeeded(reordered_y_, Y.getMap(), Y.getNumVectors(), false);
1257 Kokkos::fence();
1258
1259 for (size_t j = 0; j < X.getNumVectors(); j++) {
1260 auto X_j = X.getVector(j);
1261 auto ReorderedX_j = reordered_x_->getVectorNonConst(j);
1262 auto X_lcl = X_j->getLocalViewDevice(Tpetra::Access::ReadOnly);
1263 auto ReorderedX_lcl = ReorderedX_j->getLocalViewDevice(Tpetra::Access::ReadWrite);
1264 if (hasStreamReordered_) {
1265 local_ordinal_type stream_begin = 0;
1266 local_ordinal_type stream_end;
1267 for (int i = 0; i < num_streams_; i++) {
1268 auto perm_i = perm_v_[i];
1269 stream_end = stream_begin + perm_i.extent(0);
1270 auto X_lcl_sub = Kokkos::subview(X_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1271 auto ReorderedX_lcl_sub = Kokkos::subview(ReorderedX_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1272 Kokkos::parallel_for(
1273 Kokkos::RangePolicy<execution_space>(exec_space_instances_[i], 0, static_cast<int>(perm_i.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1274 ReorderedX_lcl_sub(perm_i(ii)) = X_lcl_sub(ii);
1275 });
1276 stream_begin = stream_end;
1277 }
1278 } else { // hasStreamsWithRCB_
1279 auto perm = perm_rcb_;
1280 Kokkos::parallel_for(
1281 Kokkos::RangePolicy<execution_space>(0, static_cast<int>(X_lcl.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1282 ReorderedX_lcl(perm(ii), 0) = X_lcl(ii, 0);
1283 });
1284 }
1285 }
1286 Kokkos::fence();
1287
1288 if (mode == Teuchos::NO_TRANS) { // Solve L (U Y) = X for Y.
1289 // Solve L Y = X for Y.
1290 L_solver_->apply(*reordered_x_, Y, mode);
1291 // Solve U Y = Y for Y.
1292 U_solver_->apply(Y, *reordered_y_, mode);
1293 } else { // Solve U^P (L^P Y) = X for Y (where P is * or T).
1294 // Solve U^P Y = X for Y.
1295 U_solver_->apply(*reordered_x_, Y, mode);
1296 // Solve L^P Y = Y for Y.
1297 L_solver_->apply(Y, *reordered_y_, mode);
1298 }
1299
1300 for (size_t j = 0; j < Y.getNumVectors(); j++) {
1301 auto Y_j = Y.getVectorNonConst(j);
1302 auto ReorderedY_j = reordered_y_->getVector(j);
1303 auto Y_lcl = Y_j->getLocalViewDevice(Tpetra::Access::ReadWrite);
1304 auto ReorderedY_lcl = ReorderedY_j->getLocalViewDevice(Tpetra::Access::ReadOnly);
1305 if (hasStreamReordered_) {
1306 local_ordinal_type stream_begin = 0;
1307 local_ordinal_type stream_end;
1308 for (int i = 0; i < num_streams_; i++) {
1309 auto perm_i = perm_v_[i];
1310 stream_end = stream_begin + perm_i.extent(0);
1311 auto Y_lcl_sub = Kokkos::subview(Y_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1312 auto ReorderedY_lcl_sub = Kokkos::subview(ReorderedY_lcl, Kokkos::make_pair(stream_begin, stream_end), 0);
1313 Kokkos::parallel_for(
1314 Kokkos::RangePolicy<execution_space>(exec_space_instances_[i], 0, static_cast<int>(perm_i.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1315 Y_lcl_sub(ii) = ReorderedY_lcl_sub(perm_i(ii));
1316 });
1317 stream_begin = stream_end;
1318 }
1319 } else { // hasStreamsWithRCB_
1320 auto perm = perm_rcb_;
1321 Kokkos::parallel_for(
1322 Kokkos::RangePolicy<execution_space>(0, static_cast<int>(Y_lcl.extent(0))), KOKKOS_LAMBDA(const int& ii) {
1323 Y_lcl(ii, 0) = ReorderedY_lcl(perm(ii), 0);
1324 });
1325 }
1326 }
1327 Kokkos::fence();
1328 } else {
1329 if (mode == Teuchos::NO_TRANS) { // Solve L (D (U Y)) = X for Y.
1330#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1331 // NOTE (Nov-15-2022):
1332 // This is a workaround for Cuda >= 11.3 (using cusparseSpSV)
1333 // since cusparseSpSV_solve() does not support in-place computation
1334 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1335
1336 // Start by solving L Y_tmp = X for Y_tmp.
1337 L_solver_->apply(X, *Y_tmp_, mode);
1338
1339 if (!this->isKokkosKernelsSpiluk_) {
1340 // Solve D Y = Y. The operation lets us do this in place in Y, so we can
1341 // write "solve D Y = Y for Y."
1342 Y_tmp_->elementWiseMultiply(one, *D_, *Y_tmp_, zero);
1343 }
1344
1345 U_solver_->apply(*Y_tmp_, Y, mode); // Solve U Y = Y_tmp.
1346#else
1347 // Start by solving L Y = X for Y.
1348 L_solver_->apply(X, Y, mode);
1349
1350 if (!this->isKokkosKernelsSpiluk_) {
1351 // Solve D Y = Y. The operation lets us do this in place in Y, so we can
1352 // write "solve D Y = Y for Y."
1353 Y.elementWiseMultiply(one, *D_, Y, zero);
1354 }
1355
1356 U_solver_->apply(Y, Y, mode); // Solve U Y = Y.
1357#endif
1358 } else { // Solve U^P (D^P (L^P Y)) = X for Y (where P is * or T).
1359#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1360 // NOTE (Nov-15-2022):
1361 // This is a workaround for Cuda >= 11.3 (using cusparseSpSV)
1362 // since cusparseSpSV_solve() does not support in-place computation
1363 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1364
1365 // Start by solving U^P Y_tmp = X for Y_tmp.
1366 U_solver_->apply(X, *Y_tmp_, mode);
1367
1368 if (!this->isKokkosKernelsSpiluk_) {
1369 // Solve D^P Y = Y.
1370 //
1371 // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we
1372 // need to do an elementwise multiply with the conjugate of
1373 // D_, not just with D_ itself.
1374 Y_tmp_->elementWiseMultiply(one, *D_, *Y_tmp_, zero);
1375 }
1376
1377 L_solver_->apply(*Y_tmp_, Y, mode); // Solve L^P Y = Y_tmp.
1378#else
1379 // Start by solving U^P Y = X for Y.
1380 U_solver_->apply(X, Y, mode);
1381
1382 if (!this->isKokkosKernelsSpiluk_) {
1383 // Solve D^P Y = Y.
1384 //
1385 // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we
1386 // need to do an elementwise multiply with the conjugate of
1387 // D_, not just with D_ itself.
1388 Y.elementWiseMultiply(one, *D_, Y, zero);
1389 }
1390
1391 L_solver_->apply(Y, Y, mode); // Solve L^P Y = Y.
1392#endif
1393 }
1394 }
1395 } else { // alpha != 1 or beta != 0
1396 if (alpha == zero) {
1397 // The special case for beta == 0 ensures that if Y contains Inf
1398 // or NaN values, we replace them with 0 (following BLAS
1399 // convention), rather than multiplying them by 0 to get NaN.
1400 if (beta == zero) {
1401 Y.putScalar(zero);
1402 } else {
1403 Y.scale(beta);
1404 }
1405 } else { // alpha != zero
1406 Impl::resetMultiVecIfNeeded(Y_tmp_, Y.getMap(), Y.getNumVectors(), false);
1407 apply(X, *Y_tmp_, mode);
1408 Y.update(alpha, *Y_tmp_, beta);
1409 }
1410 }
1411 } // end timing
1412
1413#ifdef HAVE_IFPACK2_DEBUG
1414 {
1415 Teuchos::Array<magnitude_type> norms(Y.getNumVectors());
1416 Y.norm1(norms());
1417 bool good = true;
1418 for (size_t j = 0; j < Y.getNumVectors(); ++j) {
1419 if (STM::isnaninf(norms[j])) {
1420 good = false;
1421 break;
1422 }
1423 }
1424 TEUCHOS_TEST_FOR_EXCEPTION(!good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the output Y is NaN or Inf.");
1425 }
1426#endif // HAVE_IFPACK2_DEBUG
1427
1428 ++numApply_;
1429 applyTime_ += (timer.wallTime() - startTime);
1430}
1431
1432// VINH comment out since multiply() is not needed anywhere
1433// template<class MatrixType>
1434// void RILUK<MatrixType>::
1435// multiply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1436// Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1437// const Teuchos::ETransp mode) const
1438//{
1439// const scalar_type zero = STS::zero ();
1440// const scalar_type one = STS::one ();
1441//
1442// if (mode != Teuchos::NO_TRANS) {
1443// U_->apply (X, Y, mode); //
1444// Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1445//
1446// // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we need
1447// // to do an elementwise multiply with the conjugate of D_, not
1448// // just with D_ itself.
1449// Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1450//
1451// MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y
1452// L_->apply (Y_tmp, Y, mode);
1453// Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1454// }
1455// else {
1456// L_->apply (X, Y, mode);
1457// Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1458// Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1459// MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y1
1460// U_->apply (Y_tmp, Y, mode);
1461// Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1462// }
1463// }
1464
1465template <class MatrixType>
1467 std::ostringstream os;
1468
1469 // Output is a valid YAML dictionary in flow style. If you don't
1470 // like everything on a single line, you should call describe()
1471 // instead.
1472 os << "\"Ifpack2::RILUK\": {";
1473 os << "Initialized: " << (isInitialized() ? "true" : "false") << ", "
1474 << "Computed: " << (isComputed() ? "true" : "false") << ", ";
1475
1476 os << "Level-of-fill: " << getLevelOfFill() << ", ";
1477
1478 if (isKokkosKernelsSpiluk_) os << "KK-SPILUK, ";
1479 if (isKokkosKernelsStream_) os << "KK-Stream, ";
1480
1481 if (A_.is_null()) {
1482 os << "Matrix: null";
1483 } else {
1484 os << "Global matrix dimensions: ["
1485 << A_->getGlobalNumRows() << ", " << A_->getGlobalNumCols() << "]"
1486 << ", Global nnz: " << A_->getGlobalNumEntries();
1487 }
1488
1489 if (!L_solver_.is_null()) os << ", " << L_solver_->description();
1490 if (!U_solver_.is_null()) os << ", " << U_solver_->description();
1491
1492 os << "}";
1493 return os.str();
1494}
1495
1496} // namespace Ifpack2
1497
1498#define IFPACK2_RILUK_INSTANT(S, LO, GO, N) \
1499 template class Ifpack2::RILUK<Tpetra::RowMatrix<S, LO, GO, N> >;
1500
1501#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:1199
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:1466
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:1103
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