43 #ifndef IFPACK2_CRSRILUK_DEF_HPP
44 #define IFPACK2_CRSRILUK_DEF_HPP
46 #include "Ifpack2_LocalFilter.hpp"
47 #include "Tpetra_CrsMatrix.hpp"
48 #include "Ifpack2_LocalSparseTriangularSolver.hpp"
52 template<
class MatrixType>
57 isInitialized_ (false),
62 initializeTime_ (0.0),
73 template<
class MatrixType>
78 isInitialized_ (false),
83 initializeTime_ (0.0),
94 template<
class MatrixType>
97 template<
class MatrixType>
104 template<
class MatrixType>
112 if (A.getRawPtr () != A_.getRawPtr ()) {
113 isAllocated_ =
false;
114 isInitialized_ =
false;
116 A_local_ = Teuchos::null;
117 Graph_ = Teuchos::null;
124 if (! L_solver_.is_null ()) {
127 if (! U_solver_.is_null ()) {
128 U_solver_->setMatrix (Teuchos::null);
140 template<
class MatrixType>
144 TEUCHOS_TEST_FOR_EXCEPTION(
145 L_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getL: The L factor "
146 "is null. Please call initialize() (and preferably also compute()) "
147 "before calling this method. If the input matrix has not yet been set, "
148 "you must first call setMatrix() with a nonnull input matrix before you "
149 "may call initialize() or compute().");
154 template<
class MatrixType>
155 const Tpetra::Vector<typename RILUK<MatrixType>::scalar_type,
161 TEUCHOS_TEST_FOR_EXCEPTION(
162 D_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getD: The D factor "
163 "(of diagonal entries) is null. Please call initialize() (and "
164 "preferably also compute()) before calling this method. If the input "
165 "matrix has not yet been set, you must first call setMatrix() with a "
166 "nonnull input matrix before you may call initialize() or compute().");
171 template<
class MatrixType>
175 TEUCHOS_TEST_FOR_EXCEPTION(
176 U_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getU: The U factor "
177 "is null. Please call initialize() (and preferably also compute()) "
178 "before calling this method. If the input matrix has not yet been set, "
179 "you must first call setMatrix() with a nonnull input matrix before you "
180 "may call initialize() or compute().");
185 template<
class MatrixType>
187 TEUCHOS_TEST_FOR_EXCEPTION(
188 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getNodeSmootherComplexity: "
189 "The input matrix A is null. Please call setMatrix() with a nonnull "
190 "input matrix, then call compute(), before calling this method.");
192 if(!L_.is_null() && !U_.is_null())
193 return A_->getNodeNumEntries() + L_->getNodeNumEntries() + U_->getNodeNumEntries();
200 template<
class MatrixType>
201 Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
205 TEUCHOS_TEST_FOR_EXCEPTION(
206 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getDomainMap: "
207 "The matrix is null. Please call setMatrix() with a nonnull input "
208 "before calling this method.");
211 TEUCHOS_TEST_FOR_EXCEPTION(
212 Graph_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getDomainMap: "
213 "The computed graph is null. Please call initialize() before calling "
219 template<
class MatrixType>
220 Teuchos::RCP<const Tpetra::Map<typename RILUK<MatrixType>::local_ordinal_type,
224 TEUCHOS_TEST_FOR_EXCEPTION(
225 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getRangeMap: "
226 "The matrix is null. Please call setMatrix() with a nonnull input "
227 "before calling this method.");
230 TEUCHOS_TEST_FOR_EXCEPTION(
231 Graph_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getRangeMap: "
232 "The computed graph is null. Please call initialize() before calling "
238 template<
class MatrixType>
244 if (! isAllocated_) {
253 L_ = rcp (
new crs_matrix_type (Graph_->getL_Graph ()));
254 U_ = rcp (
new crs_matrix_type (Graph_->getU_Graph ()));
255 L_->setAllToScalar (STS::zero ());
256 U_->setAllToScalar (STS::zero ());
261 D_ = rcp (
new vec_type (Graph_->getL_Graph ()->getRowMap ()));
267 template<
class MatrixType>
270 setParameters (
const Teuchos::ParameterList& params)
273 using Teuchos::Exceptions::InvalidParameterName;
274 using Teuchos::Exceptions::InvalidParameterType;
290 bool gotFillLevel =
false;
292 fillLevel = params.get<
int> (
"fact: iluk level-of-fill");
295 catch (InvalidParameterType&) {
300 catch (InvalidParameterName&) {
304 if (! gotFillLevel) {
310 catch (InvalidParameterType&) {
316 if (! gotFillLevel) {
320 fillLevel = as<int> (params.get<
magnitude_type> (
"fact: iluk level-of-fill"));
323 catch (InvalidParameterType&) {
329 if (! gotFillLevel) {
333 fillLevel = as<int> (params.get<
double> (
"fact: iluk level-of-fill"));
336 catch (InvalidParameterType& e) {
346 TEUCHOS_TEST_FOR_EXCEPTION(
349 "Ifpack2::RILUK::setParameters: We should never get here! "
350 "The method should either have read the \"fact: iluk level-of-fill\" "
351 "parameter by this point, or have thrown an exception. "
352 "Please let the Ifpack2 developers know about this bug.");
360 absThresh = params.get<
magnitude_type> (
"fact: absolute threshold");
362 catch (InvalidParameterType&) {
365 absThresh = as<magnitude_type> (params.get<
double> (
"fact: absolute threshold"));
367 catch (InvalidParameterName&) {
372 relThresh = params.get<
magnitude_type> (
"fact: relative threshold");
374 catch (InvalidParameterType&) {
377 relThresh = as<magnitude_type> (params.get<
double> (
"fact: relative threshold"));
379 catch (InvalidParameterName&) {
386 catch (InvalidParameterType&) {
389 relaxValue = as<magnitude_type> (params.get<
double> (
"fact: relax value"));
391 catch (InvalidParameterName&) {
396 L_solver_->setParameters(params);
397 U_solver_->setParameters(params);
403 LevelOfFill_ = fillLevel;
410 if (absThresh != -STM::one ()) {
411 Athresh_ = absThresh;
413 if (relThresh != -STM::one ()) {
414 Rthresh_ = relThresh;
416 if (relaxValue != -STM::one ()) {
417 RelaxValue_ = relaxValue;
422 template<
class MatrixType>
423 Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
425 return Teuchos::rcp_implicit_cast<const row_matrix_type> (A_);
429 template<
class MatrixType>
430 Teuchos::RCP<const typename RILUK<MatrixType>::crs_matrix_type>
432 return Teuchos::rcp_dynamic_cast<const crs_matrix_type> (A_,
true);
436 template<
class MatrixType>
437 Teuchos::RCP<const typename RILUK<MatrixType>::row_matrix_type>
442 using Teuchos::rcp_dynamic_cast;
443 using Teuchos::rcp_implicit_cast;
448 if (A->getRowMap ()->getComm ()->getSize () == 1 ||
449 A->getRowMap ()->isSameAs (* (A->getColMap ()))) {
456 RCP<const LocalFilter<row_matrix_type> > A_lf_r =
457 rcp_dynamic_cast<
const LocalFilter<row_matrix_type> > (A);
458 if (! A_lf_r.is_null ()) {
459 return rcp_implicit_cast<const row_matrix_type> (A_lf_r);
465 return rcp (
new LocalFilter<row_matrix_type> (A));
470 template<
class MatrixType>
475 using Teuchos::rcp_const_cast;
476 using Teuchos::rcp_dynamic_cast;
477 using Teuchos::rcp_implicit_cast;
481 const char prefix[] =
"Ifpack2::RILUK::initialize: ";
483 TEUCHOS_TEST_FOR_EXCEPTION
484 (A_.is_null (), std::runtime_error, prefix <<
"The matrix is null. Please "
485 "call setMatrix() with a nonnull input before calling this method.");
486 TEUCHOS_TEST_FOR_EXCEPTION
487 (! A_->isFillComplete (), std::runtime_error, prefix <<
"The matrix is not "
488 "fill complete. You may not invoke initialize() or compute() with this "
489 "matrix until the matrix is fill complete. If your matrix is a "
490 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
491 "range Maps, if appropriate) before calling this method.");
493 Teuchos::Time timer (
"RILUK::initialize");
495 Teuchos::TimeMonitor timeMon (timer);
504 isInitialized_ =
false;
505 isAllocated_ =
false;
507 Graph_ = Teuchos::null;
509 A_local_ = makeLocalFilter (A_);
510 TEUCHOS_TEST_FOR_EXCEPTION(
511 A_local_.is_null (), std::logic_error,
"Ifpack2::RILUK::initialize: "
512 "makeLocalFilter returned null; it failed to compute A_local. "
513 "Please report this bug to the Ifpack2 developers.");
521 RCP<const crs_matrix_type> A_local_crs =
522 rcp_dynamic_cast<const crs_matrix_type> (A_local_);
523 if (A_local_crs.is_null ()) {
528 RCP<crs_matrix_type> A_local_crs_nc =
530 A_local_->getColMap (), 0));
539 import_type
import (A_local_->getRowMap (), A_local_->getRowMap ());
540 A_local_crs_nc->doImport (*A_local_,
import, Tpetra::REPLACE);
541 A_local_crs_nc->fillComplete (A_local_->getDomainMap (), A_local_->getRangeMap ());
542 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
548 Graph_->initialize ();
550 checkOrderingConsistency (*A_local_);
551 L_solver_->setMatrix (L_);
552 L_solver_->initialize ();
553 U_solver_->setMatrix (U_);
554 U_solver_->initialize ();
561 isInitialized_ =
true;
563 initializeTime_ += timer.totalElapsedTime ();
566 template<
class MatrixType>
574 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getNodeElementList();
575 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getNodeElementList();
576 bool gidsAreConsistentlyOrdered=
true;
577 global_ordinal_type indexOfInconsistentGID=0;
578 for (global_ordinal_type i=0; i<rowGIDs.size(); ++i) {
579 if (rowGIDs[i] != colGIDs[i]) {
580 gidsAreConsistentlyOrdered=
false;
581 indexOfInconsistentGID=i;
585 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==
false, std::runtime_error,
586 "The ordering of the local GIDs in the row and column maps is not the same"
587 << std::endl <<
"at index " << indexOfInconsistentGID
588 <<
". Consistency is required, as all calculations are done with"
589 << std::endl <<
"local indexing.");
592 template<
class MatrixType>
595 initAllValues (
const row_matrix_type& A)
597 using Teuchos::ArrayRCP;
601 using Teuchos::REDUCE_SUM;
602 using Teuchos::reduceAll;
603 typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type;
605 size_t NumIn = 0, NumL = 0, NumU = 0;
606 bool DiagFound =
false;
607 size_t NumNonzeroDiags = 0;
608 size_t MaxNumEntries = A.getGlobalMaxNumRowEntries();
612 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries);
613 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
614 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
615 Teuchos::Array<scalar_type> InV(MaxNumEntries);
616 Teuchos::Array<scalar_type> LV(MaxNumEntries);
617 Teuchos::Array<scalar_type> UV(MaxNumEntries);
620 const bool ReplaceValues = L_->isStaticGraph () || L_->isLocallyIndexed ();
625 L_->setAllToScalar (STS::zero ());
626 U_->setAllToScalar (STS::zero ());
629 D_->putScalar (STS::zero ());
630 ArrayRCP<scalar_type> DV = D_->get1dViewNonConst ();
632 RCP<const map_type> rowMap = L_->getRowMap ();
642 Teuchos::ArrayView<const global_ordinal_type> nodeGIDs = rowMap->getNodeElementList();
643 for (
size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
644 local_ordinal_type local_row = myRow;
648 A.getLocalRowCopy (local_row, InI(), InV(), NumIn);
656 for (
size_t j = 0; j < NumIn; ++j) {
657 const local_ordinal_type k = InI[j];
659 if (k == local_row) {
662 DV[local_row] += Rthresh_ * InV[j] + IFPACK2_SGN(InV[j]) * Athresh_;
665 TEUCHOS_TEST_FOR_EXCEPTION(
666 true, std::runtime_error,
"Ifpack2::RILUK::initAllValues: current "
667 "GID k = " << k <<
" < 0. I'm not sure why this is an error; it is "
668 "probably an artifact of the undocumented assumptions of the "
669 "original implementation (likely copied and pasted from Ifpack). "
670 "Nevertheless, the code I found here insisted on this being an error "
671 "state, so I will throw an exception here.");
673 else if (k < local_row) {
678 else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
690 DV[local_row] = Athresh_;
695 L_->replaceLocalValues(local_row, LI(0, NumL), LV(0,NumL));
699 L_->insertLocalValues(local_row, LI(0,NumL), LV(0,NumL));
705 U_->replaceLocalValues(local_row, UI(0,NumU), UV(0,NumU));
709 U_->insertLocalValues(local_row, UI(0,NumU), UV(0,NumU));
717 isInitialized_ =
true;
721 template<
class MatrixType>
725 const char prefix[] =
"Ifpack2::RILUK::compute: ";
730 TEUCHOS_TEST_FOR_EXCEPTION
731 (A_.is_null (), std::runtime_error, prefix <<
"The matrix is null. Please "
732 "call setMatrix() with a nonnull input before calling this method.");
733 TEUCHOS_TEST_FOR_EXCEPTION
734 (! A_->isFillComplete (), std::runtime_error, prefix <<
"The matrix is not "
735 "fill complete. You may not invoke initialize() or compute() with this "
736 "matrix until the matrix is fill complete. If your matrix is a "
737 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
738 "range Maps, if appropriate) before calling this method.");
740 if (! isInitialized ()) {
744 Teuchos::Time timer (
"RILUK::compute");
747 Teuchos::TimeMonitor timeMon (timer);
753 initAllValues (*A_local_);
759 const scalar_type MaxDiagonalValue = STS::one () / MinDiagonalValue;
761 size_t NumIn, NumL, NumU;
764 const size_t MaxNumEntries =
765 L_->getNodeMaxNumRowEntries () + U_->getNodeMaxNumRowEntries () + 1;
767 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries);
768 Teuchos::Array<scalar_type> InV(MaxNumEntries);
769 size_t num_cols = U_->getColMap()->getNodeNumElements();
770 Teuchos::Array<int> colflag(num_cols);
772 Teuchos::ArrayRCP<scalar_type> DV = D_->get1dViewNonConst();
778 Teuchos::ArrayView<const local_ordinal_type> UUI;
779 Teuchos::ArrayView<const scalar_type> UUV;
780 for (
size_t j = 0; j < num_cols; ++j) {
784 for (
size_t i = 0; i < L_->getNodeNumRows (); ++i) {
789 NumIn = MaxNumEntries;
790 L_->getLocalRowCopy (local_row, InI (), InV (), NumL);
793 InI[NumL] = local_row;
795 U_->getLocalRowCopy (local_row, InI (NumL+1, MaxNumEntries-NumL-1),
796 InV (NumL+1, MaxNumEntries-NumL-1), NumU);
800 for (
size_t j = 0; j < NumIn; ++j) {
806 for (
size_t jj = 0; jj < NumL; ++jj) {
812 U_->getLocalRowView(j, UUI, UUV);
815 if (RelaxValue_ == STM::zero ()) {
816 for (
size_t k = 0; k < NumUU; ++k) {
817 const int kk = colflag[UUI[k]];
822 InV[kk] -= multiplier * UUV[k];
827 for (
size_t k = 0; k < NumUU; ++k) {
831 const int kk = colflag[UUI[k]];
833 InV[kk] -= multiplier*UUV[k];
836 diagmod -= multiplier*UUV[k];
843 L_->replaceLocalValues (local_row, InI (0, NumL), InV (0, NumL));
848 if (RelaxValue_ != STM::zero ()) {
849 DV[i] += RelaxValue_*diagmod;
852 if (STS::magnitude (DV[i]) > STS::magnitude (MaxDiagonalValue)) {
853 if (STS::real (DV[i]) < STM::zero ()) {
854 DV[i] = -MinDiagonalValue;
857 DV[i] = MinDiagonalValue;
861 DV[i] = STS::one () / DV[i];
864 for (
size_t j = 0; j < NumU; ++j) {
865 InV[NumL+1+j] *= DV[i];
870 U_->replaceLocalValues (local_row, InI (NumL+1, NumU), InV (NumL+1, NumU));
874 for (
size_t j = 0; j < NumIn; ++j) {
875 colflag[InI[j]] = -1;
886 L_->fillComplete (L_->getColMap (), A_local_->getRangeMap ());
887 U_->fillComplete (A_local_->getDomainMap (), U_->getRowMap ());
890 L_solver_->setMatrix (L_);
891 L_solver_->compute ();
892 U_solver_->setMatrix (U_);
893 U_solver_->compute ();
897 computeTime_ += timer.totalElapsedTime ();
901 template<
class MatrixType>
904 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
905 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
906 Teuchos::ETransp mode,
911 using Teuchos::rcpFromRef;
913 TEUCHOS_TEST_FOR_EXCEPTION(
914 A_.is_null (), std::runtime_error,
"Ifpack2::RILUK::apply: The matrix is "
915 "null. Please call setMatrix() with a nonnull input, then initialize() "
916 "and compute(), before calling this method.");
917 TEUCHOS_TEST_FOR_EXCEPTION(
918 ! isComputed (), std::runtime_error,
919 "Ifpack2::RILUK::apply: If you have not yet called compute(), "
920 "you must call compute() before calling this method.");
921 TEUCHOS_TEST_FOR_EXCEPTION(
922 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
923 "Ifpack2::RILUK::apply: X and Y do not have the same number of columns. "
924 "X.getNumVectors() = " << X.getNumVectors ()
925 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
926 TEUCHOS_TEST_FOR_EXCEPTION(
927 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
928 "Ifpack2::RILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
929 "complex Scalar type. Please talk to the Ifpack2 developers to get this "
930 "fixed. There is a FIXME in this file about this very issue.");
931 #ifdef HAVE_IFPACK2_DEBUG
934 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.");
935 Teuchos::Array<magnitude_type> norms (X.getNumVectors ());
938 for (
size_t j = 0; j < X.getNumVectors (); ++j) {
939 if (STM::isnaninf (norms[j])) {
944 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::RILUK::apply: The 1-norm of the input X is NaN or Inf.");
946 #endif // HAVE_IFPACK2_DEBUG
951 Teuchos::Time timer (
"RILUK::apply");
953 Teuchos::TimeMonitor timeMon (timer);
954 if (alpha == one && beta == zero) {
955 if (mode == Teuchos::NO_TRANS) {
957 L_solver_->apply (X, Y, mode);
961 Y.elementWiseMultiply (one, *D_, Y, zero);
963 U_solver_->apply (Y, Y, mode);
968 U_solver_->apply (X, Y, mode);
975 Y.elementWiseMultiply (one, *D_, Y, zero);
977 L_solver_->apply (Y, Y, mode);
991 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
992 apply (X, Y_tmp, mode);
993 Y.update (alpha, Y_tmp, beta);
998 #ifdef HAVE_IFPACK2_DEBUG
1000 Teuchos::Array<magnitude_type> norms (Y.getNumVectors ());
1003 for (
size_t j = 0; j < Y.getNumVectors (); ++j) {
1004 if (STM::isnaninf (norms[j])) {
1009 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::RILUK::apply: The 1-norm of the output Y is NaN or Inf.");
1011 #endif // HAVE_IFPACK2_DEBUG
1014 applyTime_ = timer.totalElapsedTime ();
1018 template<
class MatrixType>
1020 multiply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1021 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1022 const Teuchos::ETransp mode)
const
1024 const scalar_type zero = STS::zero ();
1025 const scalar_type one = STS::one ();
1027 if (mode != Teuchos::NO_TRANS) {
1028 U_->apply (X, Y, mode);
1029 Y.update (one, X, one);
1034 Y.elementWiseMultiply (one, *D_, Y, zero);
1036 MV Y_tmp (Y, Teuchos::Copy);
1037 L_->apply (Y_tmp, Y, mode);
1038 Y.update (one, Y_tmp, one);
1041 L_->apply (X, Y, mode);
1042 Y.update (one, X, one);
1043 Y.elementWiseMultiply (one, *D_, Y, zero);
1044 MV Y_tmp (Y, Teuchos::Copy);
1045 U_->apply (Y_tmp, Y, mode);
1046 Y.update (one, Y_tmp, one);
1051 template<
class MatrixType>
1054 std::ostringstream os;
1059 os <<
"\"Ifpack2::RILUK\": {";
1060 os <<
"Initialized: " << (isInitialized () ?
"true" :
"false") <<
", "
1061 <<
"Computed: " << (isComputed () ?
"true" :
"false") <<
", ";
1063 os <<
"Level-of-fill: " << getLevelOfFill() <<
", ";
1065 if (A_.is_null ()) {
1066 os <<
"Matrix: null";
1069 os <<
"Global matrix dimensions: ["
1070 << A_->getGlobalNumRows () <<
", " << A_->getGlobalNumCols () <<
"]"
1074 if (! L_solver_.is_null ()) os <<
", " << L_solver_->description ();
1075 if (! U_solver_.is_null ()) os <<
", " << U_solver_->description ();
1083 #define IFPACK2_RILUK_INSTANT(S,LO,GO,N) \
1084 template class Ifpack2::RILUK< Tpetra::RowMatrix<S, LO, GO, N> >;