41#ifndef IFPACK2_CRSRILUK_DEF_HPP
42#define IFPACK2_CRSRILUK_DEF_HPP
44#include "Ifpack2_LocalFilter.hpp"
45#include "Tpetra_CrsMatrix.hpp"
46#include "Teuchos_StandardParameterEntryValidators.hpp"
47#include "Ifpack2_LocalSparseTriangularSolver.hpp"
48#include "Ifpack2_Details_getParamTryingTypes.hpp"
49#include "Ifpack2_Details_getCrsMatrix.hpp"
50#include "Kokkos_Sort.hpp"
51#include "KokkosSparse_Utils.hpp"
52#include "KokkosKernels_Sorting.hpp"
63 static void loadPLTypeOption (Teuchos::Array<std::string>& type_strs, Teuchos::Array<Enum>& type_enums) {
65 type_strs[0] =
"Serial";
66 type_strs[1] =
"KSPILUK";
68 type_enums[0] = Serial;
69 type_enums[1] = KSPILUK;
74template<
class MatrixType>
75RILUK<MatrixType>::RILUK (
const Teuchos::RCP<const row_matrix_type>& Matrix_in)
80 isInitialized_ (false),
85 initializeTime_ (0.0),
97template<
class MatrixType>
98RILUK<MatrixType>::RILUK (
const Teuchos::RCP<const crs_matrix_type>& Matrix_in)
102 isAllocated_ (false),
103 isInitialized_ (false),
108 initializeTime_ (0.0),
120template<
class MatrixType>
123 if (Teuchos::nonnull (KernelHandle_))
125 KernelHandle_->destroy_spiluk_handle();
129template<
class MatrixType>
130void RILUK<MatrixType>::allocateSolvers ()
133 L_solver_->setObjectLabel(
"lower");
135 U_solver_->setObjectLabel(
"upper");
138template<
class MatrixType>
146 if (A.getRawPtr () !=
A_.getRawPtr ()) {
147 isAllocated_ =
false;
148 isInitialized_ =
false;
174template<
class MatrixType>
178 TEUCHOS_TEST_FOR_EXCEPTION(
179 L_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getL: The L factor "
180 "is null. Please call initialize() (and preferably also compute()) "
181 "before calling this method. If the input matrix has not yet been set, "
182 "you must first call setMatrix() with a nonnull input matrix before you "
183 "may call initialize() or compute().");
188template<
class MatrixType>
189const Tpetra::Vector<typename RILUK<MatrixType>::scalar_type,
195 TEUCHOS_TEST_FOR_EXCEPTION(
196 D_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getD: The D factor "
197 "(of diagonal entries) is null. Please call initialize() (and "
198 "preferably also compute()) before calling this method. If the input "
199 "matrix has not yet been set, you must first call setMatrix() with a "
200 "nonnull input matrix before you may call initialize() or compute().");
205template<
class MatrixType>
209 TEUCHOS_TEST_FOR_EXCEPTION(
210 U_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getU: The U factor "
211 "is null. Please call initialize() (and preferably also compute()) "
212 "before calling this method. If the input matrix has not yet been set, "
213 "you must first call setMatrix() with a nonnull input matrix before you "
214 "may call initialize() or compute().");
219template<
class MatrixType>
221 TEUCHOS_TEST_FOR_EXCEPTION(
222 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getNodeSmootherComplexity: "
223 "The input matrix A is null. Please call setMatrix() with a nonnull "
224 "input matrix, then call compute(), before calling this method.");
226 if(!
L_.is_null() && !
U_.is_null())
227 return A_->getLocalNumEntries() +
L_->getLocalNumEntries() +
U_->getLocalNumEntries();
233template<
class MatrixType>
234Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
238 TEUCHOS_TEST_FOR_EXCEPTION(
239 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getDomainMap: "
240 "The matrix is null. Please call setMatrix() with a nonnull input "
241 "before calling this method.");
244 TEUCHOS_TEST_FOR_EXCEPTION(
245 Graph_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getDomainMap: "
246 "The computed graph is null. Please call initialize() before calling "
248 return Graph_->getL_Graph ()->getDomainMap ();
252template<
class MatrixType>
253Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
257 TEUCHOS_TEST_FOR_EXCEPTION(
258 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getRangeMap: "
259 "The matrix is null. Please call setMatrix() with a nonnull input "
260 "before calling this method.");
263 TEUCHOS_TEST_FOR_EXCEPTION(
264 Graph_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getRangeMap: "
265 "The computed graph is null. Please call initialize() before calling "
267 return Graph_->getL_Graph ()->getRangeMap ();
271template<
class MatrixType>
272void RILUK<MatrixType>::allocate_L_and_U ()
277 if (! isAllocated_) {
286 L_ = rcp (
new crs_matrix_type (Graph_->getL_Graph ()));
287 U_ = rcp (
new crs_matrix_type (Graph_->getU_Graph ()));
288 L_->setAllToScalar (STS::zero ());
289 U_->setAllToScalar (STS::zero ());
294 D_ = rcp (
new vec_type (Graph_->getL_Graph ()->getRowMap ()));
300template<
class MatrixType>
306 using Teuchos::ParameterList;
307 using Teuchos::Array;
308 using Details::getParamTryingTypes;
309 const char prefix[] =
"Ifpack2::RILUK: ";
316 double overalloc = 2.;
327 const std::string paramName (
"fact: iluk level-of-fill");
328 getParamTryingTypes<int, int, global_ordinal_type, double, float>
329 (fillLevel, params, paramName, prefix);
334 const std::string paramName (
"fact: absolute threshold");
335 getParamTryingTypes<magnitude_type, magnitude_type, double>
336 (absThresh, params, paramName, prefix);
339 const std::string paramName (
"fact: relative threshold");
340 getParamTryingTypes<magnitude_type, magnitude_type, double>
341 (relThresh, params, paramName, prefix);
344 const std::string paramName (
"fact: relax value");
345 getParamTryingTypes<magnitude_type, magnitude_type, double>
346 (relaxValue, params, paramName, prefix);
349 const std::string paramName (
"fact: iluk overalloc");
350 getParamTryingTypes<double, double>
351 (overalloc, params, paramName, prefix);
355 Details::IlukImplType::Enum ilukimplType = Details::IlukImplType::Serial;
357 static const char typeName[] =
"fact: type";
359 if ( ! params.isType<std::string>(typeName))
break;
362 Array<std::string> ilukimplTypeStrs;
363 Array<Details::IlukImplType::Enum> ilukimplTypeEnums;
364 Details::IlukImplType::loadPLTypeOption (ilukimplTypeStrs, ilukimplTypeEnums);
365 Teuchos::StringToIntegralParameterEntryValidator<Details::IlukImplType::Enum>
366 s2i(ilukimplTypeStrs (), ilukimplTypeEnums (), typeName,
false);
368 ilukimplType = s2i.getIntegralValue(params.get<std::string>(typeName));
371 if (ilukimplType == Details::IlukImplType::KSPILUK) {
372 this->isKokkosKernelsSpiluk_ =
true;
375 this->isKokkosKernelsSpiluk_ =
false;
379 L_solver_->setParameters(params);
380 U_solver_->setParameters(params);
386 LevelOfFill_ = fillLevel;
387 Overalloc_ = overalloc;
394 if (absThresh != -STM::one ()) {
395 Athresh_ = absThresh;
397 if (relThresh != -STM::one ()) {
398 Rthresh_ = relThresh;
400 if (relaxValue != -STM::one ()) {
401 RelaxValue_ = relaxValue;
406template<
class MatrixType>
407Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
409 return Teuchos::rcp_implicit_cast<const row_matrix_type> (
A_);
413template<
class MatrixType>
414Teuchos::RCP<const typename RILUK<MatrixType>::crs_matrix_type>
416 return Teuchos::rcp_dynamic_cast<const crs_matrix_type> (
A_,
true);
420template<
class MatrixType>
421Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
422RILUK<MatrixType>::makeLocalFilter (
const Teuchos::RCP<const row_matrix_type>& A)
426 using Teuchos::rcp_dynamic_cast;
427 using Teuchos::rcp_implicit_cast;
432 if (A->getRowMap ()->getComm ()->getSize () == 1 ||
433 A->getRowMap ()->isSameAs (* (A->getColMap ()))) {
440 RCP<const LocalFilter<row_matrix_type> > A_lf_r =
441 rcp_dynamic_cast<const LocalFilter<row_matrix_type> > (A);
442 if (! A_lf_r.is_null ()) {
443 return rcp_implicit_cast<const row_matrix_type> (A_lf_r);
454template<
class MatrixType>
459 using Teuchos::rcp_const_cast;
460 using Teuchos::rcp_dynamic_cast;
461 using Teuchos::rcp_implicit_cast;
462 using Teuchos::Array;
463 using Teuchos::ArrayView;
467 const char prefix[] =
"Ifpack2::RILUK::initialize: ";
469 TEUCHOS_TEST_FOR_EXCEPTION
470 (
A_.is_null (), std::runtime_error, prefix <<
"The matrix is null. Please "
471 "call setMatrix() with a nonnull input before calling this method.");
472 TEUCHOS_TEST_FOR_EXCEPTION
473 (!
A_->isFillComplete (), std::runtime_error, prefix <<
"The matrix is not "
474 "fill complete. You may not invoke initialize() or compute() with this "
475 "matrix until the matrix is fill complete. If your matrix is a "
476 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
477 "range Maps, if appropriate) before calling this method.");
479 Teuchos::Time timer (
"RILUK::initialize");
480 double startTime = timer.wallTime();
482 Teuchos::TimeMonitor timeMon (timer);
491 isInitialized_ =
false;
492 isAllocated_ =
false;
497 TEUCHOS_TEST_FOR_EXCEPTION(
498 A_local_.is_null (), std::logic_error,
"Ifpack2::RILUK::initialize: "
499 "makeLocalFilter returned null; it failed to compute A_local. "
500 "Please report this bug to the Ifpack2 developers.");
509 this->KernelHandle_ = Teuchos::rcp (
new kk_handle_type ());
510 KernelHandle_->create_spiluk_handle( KokkosSparse::Experimental::SPILUKAlgorithm::SEQLVLSCHD_TP1,
512 2*
A_local_->getLocalNumEntries()*(LevelOfFill_+1),
513 2*
A_local_->getLocalNumEntries()*(LevelOfFill_+1) );
517 RCP<const crs_matrix_type> A_local_crs = Details::getCrsMatrix(A_local_);
518 if(A_local_crs.is_null()) {
520 Array<size_t> entriesPerRow(numRows);
522 entriesPerRow[i] =
A_local_->getNumEntriesInLocalRow(i);
524 RCP<crs_matrix_type> A_local_crs_nc =
525 rcp (
new crs_matrix_type (A_local_->getRowMap (),
526 A_local_->getColMap (),
529 nonconst_local_inds_host_view_type indices(
"indices",A_local_->getLocalMaxNumRowEntries());
530 nonconst_values_host_view_type values(
"values",A_local_->getLocalMaxNumRowEntries());
532 size_t numEntries = 0;
533 A_local_->getLocalRowCopy(i, indices, values, numEntries);
534 A_local_crs_nc->insertLocalValues(i, numEntries,
reinterpret_cast<scalar_type*
>(values.data()), indices.data());
536 A_local_crs_nc->fillComplete (A_local_->getDomainMap (), A_local_->getRangeMap ());
537 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
539 Graph_ = rcp (
new Ifpack2::IlukGraph<crs_graph_type,kk_handle_type> (A_local_crs->getCrsGraph (),
540 LevelOfFill_, 0, Overalloc_));
544 if (this->isKokkosKernelsSpiluk_) Graph_->initialize (KernelHandle_);
545 else Graph_->initialize ();
548 checkOrderingConsistency (*A_local_);
549 L_solver_->setMatrix (L_);
550 L_solver_->initialize ();
554#if !defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) || !defined(KOKKOS_ENABLE_CUDA) || (CUDA_VERSION < 11030)
555 L_solver_->compute ();
557 U_solver_->setMatrix (U_);
559#if !defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) || !defined(KOKKOS_ENABLE_CUDA) || (CUDA_VERSION < 11030)
568 isInitialized_ =
true;
570 initializeTime_ += (timer.wallTime() - startTime);
573template<
class MatrixType>
581 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getLocalElementList();
582 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getLocalElementList();
583 bool gidsAreConsistentlyOrdered=
true;
586 if (rowGIDs[i] != colGIDs[i]) {
587 gidsAreConsistentlyOrdered=
false;
588 indexOfInconsistentGID=i;
592 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==
false, std::runtime_error,
593 "The ordering of the local GIDs in the row and column maps is not the same"
594 << std::endl <<
"at index " << indexOfInconsistentGID
595 <<
". Consistency is required, as all calculations are done with"
596 << std::endl <<
"local indexing.");
599template<
class MatrixType>
604 using Teuchos::ArrayRCP;
608 using Teuchos::REDUCE_SUM;
609 using Teuchos::reduceAll;
610 typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type>
map_type;
612 size_t NumIn = 0, NumL = 0, NumU = 0;
613 bool DiagFound =
false;
614 size_t NumNonzeroDiags = 0;
615 size_t MaxNumEntries = A.getGlobalMaxNumRowEntries();
619 nonconst_local_inds_host_view_type InI(
"InI",MaxNumEntries);
620 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
621 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
622 nonconst_values_host_view_type InV(
"InV",MaxNumEntries);
623 Teuchos::Array<scalar_type> LV(MaxNumEntries);
624 Teuchos::Array<scalar_type> UV(MaxNumEntries);
627 const bool ReplaceValues = L_->isStaticGraph () || L_->isLocallyIndexed ();
632 L_->setAllToScalar (STS::zero ());
633 U_->setAllToScalar (STS::zero ());
636 D_->putScalar (STS::zero ());
637 auto DV = Kokkos::subview(D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
639 RCP<const map_type> rowMap = L_->getRowMap ();
649 Teuchos::ArrayView<const global_ordinal_type> nodeGIDs = rowMap->getLocalElementList();
650 for (
size_t myRow=0; myRow<A.getLocalNumRows(); ++myRow) {
655 A.getLocalRowCopy (local_row, InI, InV, NumIn);
663 for (
size_t j = 0; j < NumIn; ++j) {
666 if (k == local_row) {
669 DV(local_row) += Rthresh_ * InV[j] + IFPACK2_SGN(InV[j]) * Athresh_;
672 TEUCHOS_TEST_FOR_EXCEPTION(
673 true, std::runtime_error,
"Ifpack2::RILUK::initAllValues: current "
674 "GID k = " << k <<
" < 0. I'm not sure why this is an error; it is "
675 "probably an artifact of the undocumented assumptions of the "
676 "original implementation (likely copied and pasted from Ifpack). "
677 "Nevertheless, the code I found here insisted on this being an error "
678 "state, so I will throw an exception here.");
680 else if (k < local_row) {
685 else if (Teuchos::as<size_t>(k) <= rowMap->getLocalNumElements()) {
697 DV(local_row) = Athresh_;
702 L_->replaceLocalValues(local_row, LI(0, NumL), LV(0,NumL));
706 L_->insertLocalValues(local_row, LI(0,NumL), LV(0,NumL));
712 U_->replaceLocalValues(local_row, UI(0,NumU), UV(0,NumU));
716 U_->insertLocalValues(local_row, UI(0,NumU), UV(0,NumU));
724 isInitialized_ =
true;
728template<
class MatrixType>
733 using Teuchos::rcp_const_cast;
734 using Teuchos::rcp_dynamic_cast;
735 using Teuchos::Array;
736 using Teuchos::ArrayView;
737 const char prefix[] =
"Ifpack2::RILUK::compute: ";
742 TEUCHOS_TEST_FOR_EXCEPTION
743 (
A_.is_null (), std::runtime_error, prefix <<
"The matrix is null. Please "
744 "call setMatrix() with a nonnull input before calling this method.");
745 TEUCHOS_TEST_FOR_EXCEPTION
746 (!
A_->isFillComplete (), std::runtime_error, prefix <<
"The matrix is not "
747 "fill complete. You may not invoke initialize() or compute() with this "
748 "matrix until the matrix is fill complete. If your matrix is a "
749 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
750 "range Maps, if appropriate) before calling this method.");
756 Teuchos::Time timer (
"RILUK::compute");
759 Teuchos::TimeMonitor timeMon (timer);
760 double startTime = timer.wallTime();
773 const scalar_type MaxDiagonalValue = STS::one () / MinDiagonalValue;
775 size_t NumIn, NumL, NumU;
778 const size_t MaxNumEntries =
779 L_->getLocalMaxNumRowEntries () +
U_->getLocalMaxNumRowEntries () + 1;
781 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries);
782 Teuchos::Array<scalar_type> InV(MaxNumEntries);
783 size_t num_cols =
U_->getColMap()->getLocalNumElements();
784 Teuchos::Array<int> colflag(num_cols);
786 auto DV = Kokkos::subview(
D_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
790 for (
size_t j = 0; j < num_cols; ++j) {
793 using IST =
typename row_matrix_type::impl_scalar_type;
794 for (
size_t i = 0; i <
L_->getLocalNumRows (); ++i) {
798 local_inds_host_view_type UUI;
799 values_host_view_type UUV;
803 NumIn = MaxNumEntries;
804 nonconst_local_inds_host_view_type InI_v(InI.data(),MaxNumEntries);
805 nonconst_values_host_view_type InV_v(
reinterpret_cast<IST*
>(InV.data()),MaxNumEntries);
807 L_->getLocalRowCopy (local_row, InI_v , InV_v, NumL);
810 InI[NumL] = local_row;
812 nonconst_local_inds_host_view_type InI_sub(InI.data()+NumL+1,MaxNumEntries-NumL-1);
813 nonconst_values_host_view_type InV_sub(
reinterpret_cast<IST*
>(InV.data())+NumL+1,MaxNumEntries-NumL-1);
815 U_->getLocalRowCopy (local_row, InI_sub,InV_sub, NumU);
819 for (
size_t j = 0; j < NumIn; ++j) {
825 for (
size_t jj = 0; jj < NumL; ++jj) {
827 IST multiplier = InV[jj];
831 U_->getLocalRowView(j, UUI, UUV);
834 if (RelaxValue_ == STM::zero ()) {
835 for (
size_t k = 0; k < NumUU; ++k) {
836 const int kk = colflag[UUI[k]];
841 InV[kk] -=
static_cast<scalar_type>(multiplier * UUV[k]);
847 for (
size_t k = 0; k < NumUU; ++k) {
851 const int kk = colflag[UUI[k]];
853 InV[kk] -=
static_cast<scalar_type>(multiplier*UUV[k]);
856 diagmod -=
static_cast<scalar_type>(multiplier*UUV[k]);
864 L_->replaceLocalValues (local_row, InI (0, NumL), InV (0, NumL));
869 if (RelaxValue_ != STM::zero ()) {
870 DV(i) += RelaxValue_*diagmod;
873 if (STS::magnitude (DV(i)) > STS::magnitude (MaxDiagonalValue)) {
874 if (STS::real (DV(i)) < STM::zero ()) {
875 DV(i) = -MinDiagonalValue;
878 DV(i) = MinDiagonalValue;
885 for (
size_t j = 0; j < NumU; ++j) {
891 U_->replaceLocalValues (local_row, InI (NumL+1, NumU), InV (NumL+1, NumU));
895 for (
size_t j = 0; j < NumIn; ++j) {
896 colflag[InI[j]] = -1;
907 L_->fillComplete (
L_->getColMap (),
A_local_->getRangeMap ());
908 U_->fillComplete (
A_local_->getDomainMap (),
U_->getRowMap ());
918 RCP<const crs_matrix_type> A_local_crs = Details::getCrsMatrix(
A_local_);
919 if(A_local_crs.is_null()) {
921 Array<size_t> entriesPerRow(numRows);
923 entriesPerRow[i] =
A_local_->getNumEntriesInLocalRow(i);
925 RCP<crs_matrix_type> A_local_crs_nc =
930 nonconst_local_inds_host_view_type indices(
"indices",
A_local_->getLocalMaxNumRowEntries());
931 nonconst_values_host_view_type values(
"values",
A_local_->getLocalMaxNumRowEntries());
933 size_t numEntries = 0;
934 A_local_->getLocalRowCopy(i, indices, values, numEntries);
935 A_local_crs_nc->insertLocalValues(i, numEntries,
reinterpret_cast<scalar_type*
>(values.data()),indices.data());
937 A_local_crs_nc->fillComplete (
A_local_->getDomainMap (),
A_local_->getRangeMap ());
938 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
940 auto lclMtx = A_local_crs->getLocalMatrixDevice();
941 A_local_rowmap_ = lclMtx.graph.row_map;
942 A_local_entries_ = lclMtx.graph.entries;
943 A_local_values_ = lclMtx.values;
949 if (
L_->isStaticGraph () ||
L_->isLocallyIndexed ()) {
950 L_->setAllToScalar (STS::zero ());
951 U_->setAllToScalar (STS::zero ());
954 using row_map_type =
typename crs_matrix_type::local_matrix_device_type::row_map_type;
957 auto lclL =
L_->getLocalMatrixDevice();
958 row_map_type L_rowmap = lclL.graph.row_map;
959 auto L_entries = lclL.graph.entries;
960 auto L_values = lclL.values;
962 auto lclU =
U_->getLocalMatrixDevice();
963 row_map_type U_rowmap = lclU.graph.row_map;
964 auto U_entries = lclU.graph.entries;
965 auto U_values = lclU.values;
967 KokkosSparse::Experimental::spiluk_numeric( KernelHandle_.getRawPtr(), LevelOfFill_,
968 A_local_rowmap_, A_local_entries_, A_local_values_,
969 L_rowmap, L_entries, L_values, U_rowmap, U_entries, U_values );
972 L_->fillComplete (
L_->getColMap (),
A_local_->getRangeMap ());
973 U_->fillComplete (
A_local_->getDomainMap (),
U_->getRowMap ());
983 computeTime_ += (timer.wallTime() - startTime);
987template<
class MatrixType>
990apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
991 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
992 Teuchos::ETransp mode,
997 using Teuchos::rcpFromRef;
999 TEUCHOS_TEST_FOR_EXCEPTION(
1000 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::apply: The matrix is "
1001 "null. Please call setMatrix() with a nonnull input, then initialize() "
1002 "and compute(), before calling this method.");
1003 TEUCHOS_TEST_FOR_EXCEPTION(
1005 "Ifpack2::RILUK::apply: If you have not yet called compute(), "
1006 "you must call compute() before calling this method.");
1007 TEUCHOS_TEST_FOR_EXCEPTION(
1008 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
1009 "Ifpack2::RILUK::apply: X and Y do not have the same number of columns. "
1010 "X.getNumVectors() = " << X.getNumVectors ()
1011 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
1012 TEUCHOS_TEST_FOR_EXCEPTION(
1013 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
1014 "Ifpack2::RILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
1015 "complex Scalar type. Please talk to the Ifpack2 developers to get this "
1016 "fixed. There is a FIXME in this file about this very issue.");
1017#ifdef HAVE_IFPACK2_DEBUG
1020 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.");
1021 Teuchos::Array<magnitude_type> norms (X.getNumVectors ());
1024 for (
size_t j = 0; j < X.getNumVectors (); ++j) {
1025 if (STM::isnaninf (norms[j])) {
1030 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::RILUK::apply: The 1-norm of the input X is NaN or Inf.");
1037 Teuchos::Time timer (
"RILUK::apply");
1038 double startTime = timer.wallTime();
1040 Teuchos::TimeMonitor timeMon (timer);
1041 if (alpha == one && beta == zero) {
1042 if (mode == Teuchos::NO_TRANS) {
1043#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1047 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
1055 Y_tmp.elementWiseMultiply (one, *
D_, Y_tmp, zero);
1066 Y.elementWiseMultiply (one, *
D_, Y, zero);
1073#if defined(KOKKOSKERNELS_ENABLE_TPL_CUSPARSE) && defined(KOKKOS_ENABLE_CUDA) && (CUDA_VERSION >= 11030)
1077 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
1088 Y_tmp.elementWiseMultiply (one, *
D_, Y_tmp, zero);
1102 Y.elementWiseMultiply (one, *
D_, Y, zero);
1110 if (alpha == zero) {
1120 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
1121 apply (X, Y_tmp, mode);
1122 Y.update (alpha, Y_tmp, beta);
1127#ifdef HAVE_IFPACK2_DEBUG
1129 Teuchos::Array<magnitude_type> norms (Y.getNumVectors ());
1132 for (
size_t j = 0; j < Y.getNumVectors (); ++j) {
1133 if (STM::isnaninf (norms[j])) {
1138 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::RILUK::apply: The 1-norm of the output Y is NaN or Inf.");
1143 applyTime_ += (timer.wallTime() - startTime);
1180template<
class MatrixType>
1183 std::ostringstream os;
1188 os <<
"\"Ifpack2::RILUK\": {";
1189 os <<
"Initialized: " << (
isInitialized () ?
"true" :
"false") <<
", "
1190 <<
"Computed: " << (
isComputed () ?
"true" :
"false") <<
", ";
1194 if (
A_.is_null ()) {
1195 os <<
"Matrix: null";
1198 os <<
"Global matrix dimensions: ["
1199 <<
A_->getGlobalNumRows () <<
", " <<
A_->getGlobalNumCols () <<
"]"
1200 <<
", Global nnz: " <<
A_->getGlobalNumEntries();
1212#define IFPACK2_RILUK_INSTANT(S,LO,GO,N) \
1213 template class Ifpack2::RILUK< Tpetra::RowMatrix<S, LO, GO, N> >;
Access only local rows and columns of a sparse matrix.
Definition Ifpack2_LocalFilter_decl.hpp:163
"Preconditioner" that solves local sparse triangular systems.
Definition Ifpack2_LocalSparseTriangularSolver_decl.hpp:88
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition Ifpack2_RILUK_decl.hpp:254
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector).
Definition Ifpack2_RILUK_decl.hpp:293
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:455
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:290
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition Ifpack2_RILUK_def.hpp:121
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:990
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_RILUK_decl.hpp:269
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:256
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition Ifpack2_RILUK_decl.hpp:592
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:611
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:207
std::string description() const
A one-line description of this object.
Definition Ifpack2_RILUK_def.hpp:1181
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_RILUK_def.hpp:140
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:257
bool isComputed() const
Whether compute() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:373
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:193
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:266
Teuchos::RCP< const row_matrix_type > A_local_
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_t...
Definition Ifpack2_RILUK_decl.hpp:601
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_RILUK_def.hpp:408
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition Ifpack2_RILUK_decl.hpp:609
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type >, kk_handle_type > > Graph_
The ILU(k) graph.
Definition Ifpack2_RILUK_decl.hpp:597
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition Ifpack2_RILUK_decl.hpp:535
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:415
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:237
void setParameters(const Teuchos::ParameterList ¶ms)
Definition Ifpack2_RILUK_def.hpp:303
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:260
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:607
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:263
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:176
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition Ifpack2_RILUK_decl.hpp:613
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_RILUK_def.hpp:220
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:729
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition Ifpack2_RILUK_decl.hpp:616
bool isKokkosKernelsSpiluk_
Optional KokkosKernels implementation.
Definition Ifpack2_RILUK_decl.hpp:638
bool isInitialized() const
Whether initialize() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:369
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74