41#ifndef IFPACK2_RELAXATION_DEF_HPP
42#define IFPACK2_RELAXATION_DEF_HPP
44#include "Teuchos_StandardParameterEntryValidators.hpp"
45#include "Teuchos_TimeMonitor.hpp"
46#include "Tpetra_CrsMatrix.hpp"
47#include "Tpetra_BlockCrsMatrix.hpp"
48#include "Tpetra_BlockView.hpp"
50#include "Ifpack2_Details_getCrsMatrix.hpp"
51#include "MatrixMarket_Tpetra.hpp"
52#include "Tpetra_Details_residual.hpp"
56#include "KokkosSparse_gauss_seidel.hpp"
60 class NonnegativeIntValidator :
public Teuchos::ParameterEntryValidator {
63 NonnegativeIntValidator () {}
66 Teuchos::ParameterEntryValidator::ValidStringsList validStringValues ()
const {
72 validate (
const Teuchos::ParameterEntry& entry,
73 const std::string& paramName,
74 const std::string& sublistName)
const
77 Teuchos::any anyVal = entry.getAny (
true);
78 const std::string entryName = entry.getAny (
false).typeName ();
80 TEUCHOS_TEST_FOR_EXCEPTION(
81 anyVal.type () != typeid (
int),
82 Teuchos::Exceptions::InvalidParameterType,
83 "Parameter \"" << paramName <<
"\" in sublist \"" << sublistName
84 <<
"\" has the wrong type." << endl <<
"Parameter: " << paramName
85 << endl <<
"Type specified: " << entryName << endl
86 <<
"Type required: int" << endl);
88 const int val = Teuchos::any_cast<int> (anyVal);
89 TEUCHOS_TEST_FOR_EXCEPTION(
90 val < 0, Teuchos::Exceptions::InvalidParameterValue,
91 "Parameter \"" << paramName <<
"\" in sublist \"" << sublistName
92 <<
"\" is negative." << endl <<
"Parameter: " << paramName
93 << endl <<
"Value specified: " << val << endl
94 <<
"Required range: [0, INT_MAX]" << endl);
98 const std::string getXMLTypeName ()
const {
99 return "NonnegativeIntValidator";
104 printDoc (
const std::string& docString,
105 std::ostream &out)
const
107 Teuchos::StrUtils::printLines (out,
"# ", docString);
108 out <<
"#\tValidator Used: " << std::endl;
109 out <<
"#\t\tNonnegativeIntValidator" << std::endl;
116 template<class Scalar, const bool isOrdinal=Teuchos::ScalarTraits<Scalar>::isOrdinal>
120 static const Scalar eps ();
124 template<
class Scalar>
125 class SmallTraits<Scalar, true> {
127 static const Scalar eps () {
128 return Teuchos::ScalarTraits<Scalar>::one ();
133 template<
class Scalar>
134 class SmallTraits<Scalar, false> {
136 static const Scalar eps () {
137 return Teuchos::ScalarTraits<Scalar>::eps ();
142 template<
class Scalar,
143 const bool isComplex = Teuchos::ScalarTraits<Scalar>::isComplex>
144 struct RealTraits {};
146 template<
class Scalar>
147 struct RealTraits<Scalar, false> {
148 using val_type = Scalar;
149 using mag_type = Scalar;
150 static KOKKOS_INLINE_FUNCTION mag_type real (
const val_type& z) {
155 template<
class Scalar>
156 struct RealTraits<Scalar, true> {
157 using val_type = Scalar;
158 using mag_type =
typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
159 static KOKKOS_INLINE_FUNCTION mag_type real (
const val_type& z) {
160 return Kokkos::ArithTraits<val_type>::real (z);
164 template<
class Scalar>
165 KOKKOS_INLINE_FUNCTION
typename RealTraits<Scalar>::mag_type
166 getRealValue (
const Scalar& z) {
167 return RealTraits<Scalar>::real (z);
174template<
class MatrixType>
177updateCachedMultiVector (
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type>>& map,
178 size_t numVecs)
const
182 if (cachedMV_.is_null () ||
183 map.get () != cachedMV_->getMap ().get () ||
184 cachedMV_->getNumVectors () != numVecs) {
187 cachedMV_ = Teuchos::rcp (
new MV (map, numVecs,
false));
191template<
class MatrixType>
193setMatrix(
const Teuchos::RCP<const row_matrix_type>& A)
195 if (A.getRawPtr() != A_.getRawPtr()) {
196 Importer_ = Teuchos::null;
197 pointImporter_ = Teuchos::null;
198 Diagonal_ = Teuchos::null;
199 isInitialized_ =
false;
201 diagOffsets_ = Kokkos::View<size_t*, device_type>();
202 savedDiagOffsets_ =
false;
203 hasBlockCrsMatrix_ =
false;
204 serialGaussSeidel_ = Teuchos::null;
205 if (! A.is_null ()) {
206 IsParallel_ = (A->getRowMap ()->
getComm ()->getSize () > 1);
212template<
class MatrixType>
214Relaxation (
const Teuchos::RCP<const row_matrix_type>& A)
216 IsParallel_ ((A.is_null () || A->getRowMap ().is_null () || A->getRowMap ()->
getComm ().is_null ()) ?
218 A->getRowMap ()->
getComm ()->getSize () > 1)
220 this->setObjectLabel (
"Ifpack2::Relaxation");
224template<
class MatrixType>
225Teuchos::RCP<const Teuchos::ParameterList>
228 using Teuchos::Array;
229 using Teuchos::ParameterList;
230 using Teuchos::parameterList;
233 using Teuchos::rcp_const_cast;
234 using Teuchos::rcp_implicit_cast;
235 using Teuchos::setStringToIntegralParameter;
236 typedef Teuchos::ParameterEntryValidator PEV;
238 if (validParams_.is_null ()) {
239 RCP<ParameterList> pl = parameterList (
"Ifpack2::Relaxation");
243 Array<std::string> precTypes (8);
244 precTypes[0] =
"Jacobi";
245 precTypes[1] =
"Gauss-Seidel";
246 precTypes[2] =
"Symmetric Gauss-Seidel";
247 precTypes[3] =
"MT Gauss-Seidel";
248 precTypes[4] =
"MT Symmetric Gauss-Seidel";
249 precTypes[5] =
"Richardson";
250 precTypes[6] =
"Two-stage Gauss-Seidel";
251 precTypes[7] =
"Two-stage Symmetric Gauss-Seidel";
252 Array<Details::RelaxationType> precTypeEnums (8);
253 precTypeEnums[0] = Details::JACOBI;
254 precTypeEnums[1] = Details::GS;
255 precTypeEnums[2] = Details::SGS;
256 precTypeEnums[3] = Details::MTGS;
257 precTypeEnums[4] = Details::MTSGS;
258 precTypeEnums[5] = Details::RICHARDSON;
259 precTypeEnums[6] = Details::GS2;
260 precTypeEnums[7] = Details::SGS2;
261 const std::string defaultPrecType (
"Jacobi");
262 setStringToIntegralParameter<Details::RelaxationType> (
"relaxation: type",
263 defaultPrecType,
"Relaxation method", precTypes (), precTypeEnums (),
266 const int numSweeps = 1;
267 RCP<PEV> numSweepsValidator =
268 rcp_implicit_cast<PEV> (rcp (
new NonnegativeIntValidator));
269 pl->set (
"relaxation: sweeps", numSweeps,
"Number of relaxation sweeps",
270 rcp_const_cast<const PEV> (numSweepsValidator));
273 const int numOuterSweeps = 1;
274 RCP<PEV> numOuterSweepsValidator =
275 rcp_implicit_cast<PEV> (rcp (
new NonnegativeIntValidator));
276 pl->set (
"relaxation: outer sweeps", numOuterSweeps,
"Number of outer local relaxation sweeps for two-stage GS",
277 rcp_const_cast<const PEV> (numOuterSweepsValidator));
279 const int numInnerSweeps = 1;
280 RCP<PEV> numInnerSweepsValidator =
281 rcp_implicit_cast<PEV> (rcp (
new NonnegativeIntValidator));
282 pl->set (
"relaxation: inner sweeps", numInnerSweeps,
"Number of inner local relaxation sweeps for two-stage GS",
283 rcp_const_cast<const PEV> (numInnerSweepsValidator));
285 const scalar_type innerDampingFactor = STS::one ();
286 pl->set (
"relaxation: inner damping factor", innerDampingFactor,
"Damping factor for the inner sweep of two-stage GS");
288 const bool innerSpTrsv =
false;
289 pl->set (
"relaxation: inner sparse-triangular solve", innerSpTrsv,
"Specify whether to use sptrsv instead of JR iterations for two-stage GS");
291 const bool compactForm =
false;
292 pl->set (
"relaxation: compact form", compactForm,
"Specify whether to use compact form of recurrence for two-stage GS");
295 pl->set (
"relaxation: damping factor", dampingFactor);
297 const bool zeroStartingSolution =
true;
298 pl->set (
"relaxation: zero starting solution", zeroStartingSolution);
300 const bool doBackwardGS =
false;
301 pl->set (
"relaxation: backward mode", doBackwardGS);
303 const bool doL1Method =
false;
304 pl->set (
"relaxation: use l1", doL1Method);
306 const magnitude_type l1eta = (STM::one() + STM::one() + STM::one()) /
307 (STM::one() + STM::one());
308 pl->set (
"relaxation: l1 eta", l1eta);
311 pl->set (
"relaxation: min diagonal value", minDiagonalValue);
313 const bool fixTinyDiagEntries =
false;
314 pl->set (
"relaxation: fix tiny diagonal entries", fixTinyDiagEntries);
316 const bool checkDiagEntries =
false;
317 pl->set (
"relaxation: check diagonal entries", checkDiagEntries);
319 Teuchos::ArrayRCP<local_ordinal_type> localSmoothingIndices = Teuchos::null;
320 pl->set(
"relaxation: local smoothing indices", localSmoothingIndices);
322 const bool is_matrix_structurally_symmetric =
false;
323 pl->set(
"relaxation: symmetric matrix structure", is_matrix_structurally_symmetric);
325 const bool ifpack2_dump_matrix =
false;
326 pl->set(
"relaxation: ifpack2 dump matrix", ifpack2_dump_matrix);
328 const int cluster_size = 1;
329 pl->set(
"relaxation: mtgs cluster size", cluster_size);
331 pl->set(
"relaxation: mtgs coloring algorithm",
"Default");
333 const int long_row_threshold = 0;
334 pl->set(
"relaxation: long row threshold", long_row_threshold);
336 const bool timer_for_apply =
true;
337 pl->set(
"timer for apply", timer_for_apply);
339 validParams_ = rcp_const_cast<const ParameterList> (pl);
345template<
class MatrixType>
346void Relaxation<MatrixType>::setParametersImpl (Teuchos::ParameterList& pl)
348 using Teuchos::getIntegralValue;
349 using Teuchos::ParameterList;
351 typedef scalar_type ST;
353 if (pl.isType<
double>(
"relaxation: damping factor")) {
356 ST df = pl.get<
double>(
"relaxation: damping factor");
357 pl.remove(
"relaxation: damping factor");
358 pl.set(
"relaxation: damping factor",df);
363 const Details::RelaxationType precType =
364 getIntegralValue<Details::RelaxationType> (pl,
"relaxation: type");
365 const int numSweeps = pl.get<
int> (
"relaxation: sweeps");
366 const ST dampingFactor = pl.get<ST> (
"relaxation: damping factor");
367 const bool zeroStartSol = pl.get<
bool> (
"relaxation: zero starting solution");
368 const bool doBackwardGS = pl.get<
bool> (
"relaxation: backward mode");
369 const bool doL1Method = pl.get<
bool> (
"relaxation: use l1");
370 const magnitude_type l1Eta = pl.get<magnitude_type> (
"relaxation: l1 eta");
371 const ST minDiagonalValue = pl.get<ST> (
"relaxation: min diagonal value");
372 const bool fixTinyDiagEntries = pl.get<
bool> (
"relaxation: fix tiny diagonal entries");
373 const bool checkDiagEntries = pl.get<
bool> (
"relaxation: check diagonal entries");
374 const bool is_matrix_structurally_symmetric = pl.get<
bool> (
"relaxation: symmetric matrix structure");
375 const bool ifpack2_dump_matrix = pl.get<
bool> (
"relaxation: ifpack2 dump matrix");
376 const bool timer_for_apply = pl.get<
bool> (
"timer for apply");
377 int cluster_size = 1;
378 if(pl.isParameter (
"relaxation: mtgs cluster size"))
379 cluster_size = pl.get<
int> (
"relaxation: mtgs cluster size");
380 int long_row_threshold = 0;
381 if(pl.isParameter (
"relaxation: long row threshold"))
382 long_row_threshold = pl.get<
int> (
"relaxation: long row threshold");
383 std::string color_algo_name = pl.get<std::string>(
"relaxation: mtgs coloring algorithm");
385 for(
char& c : color_algo_name)
387 if(color_algo_name ==
"default")
388 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_DEFAULT;
389 else if(color_algo_name ==
"serial")
390 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_SERIAL;
391 else if(color_algo_name ==
"vb")
392 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_VB;
393 else if(color_algo_name ==
"vbbit")
394 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_VBBIT;
395 else if(color_algo_name ==
"vbcs")
396 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_VBCS;
397 else if(color_algo_name ==
"vbd")
398 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_VBD;
399 else if(color_algo_name ==
"vbdbit")
400 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_VBDBIT;
401 else if(color_algo_name ==
"eb")
402 this->mtColoringAlgorithm_ = KokkosGraph::COLORING_EB;
405 std::ostringstream msg;
406 msg <<
"Ifpack2::Relaxation: 'relaxation: mtgs coloring algorithm' = '" << color_algo_name <<
"' is not valid.\n";
407 msg <<
"Choices (not case sensitive) are: Default, Serial, VB, VBBIT, VBCS, VBD, VBDBIT, EB.";
408 TEUCHOS_TEST_FOR_EXCEPTION(
409 true, std::invalid_argument, msg.str());
412 Teuchos::ArrayRCP<local_ordinal_type> localSmoothingIndices = pl.get<Teuchos::ArrayRCP<local_ordinal_type> >(
"relaxation: local smoothing indices");
415 if (!std::is_same<double, ST>::value && pl.isType<
double>(
"relaxation: inner damping factor")) {
418 ST df = pl.get<
double>(
"relaxation: inner damping factor");
419 pl.remove(
"relaxation: inner damping factor");
420 pl.set(
"relaxation: inner damping factor",df);
423 if (long_row_threshold > 0) {
424 TEUCHOS_TEST_FOR_EXCEPTION(
425 cluster_size != 1, std::invalid_argument,
"Ifpack2::Relaxation: "
426 "Requested long row MTGS/MTSGS algorithm and cluster GS/SGS, but those are not compatible.");
427 TEUCHOS_TEST_FOR_EXCEPTION(
428 precType != Details::RelaxationType::MTGS && precType != Details::RelaxationType::MTSGS,
429 std::invalid_argument,
"Ifpack2::Relaxation: "
430 "Requested long row MTGS/MTSGS algorithm, but this is only compatible with preconditioner types "
431 "'MT Gauss-Seidel' and 'MT Symmetric Gauss-Seidel'.");
434 const ST innerDampingFactor = pl.get<ST> (
"relaxation: inner damping factor");
435 const int numInnerSweeps = pl.get<
int> (
"relaxation: inner sweeps");
436 const int numOuterSweeps = pl.get<
int> (
"relaxation: outer sweeps");
437 const bool innerSpTrsv = pl.get<
bool> (
"relaxation: inner sparse-triangular solve");
438 const bool compactForm = pl.get<
bool> (
"relaxation: compact form");
441 PrecType_ = precType;
442 NumSweeps_ = numSweeps;
443 DampingFactor_ = dampingFactor;
444 ZeroStartingSolution_ = zeroStartSol;
445 DoBackwardGS_ = doBackwardGS;
446 DoL1Method_ = doL1Method;
448 MinDiagonalValue_ = minDiagonalValue;
449 fixTinyDiagEntries_ = fixTinyDiagEntries;
450 checkDiagEntries_ = checkDiagEntries;
451 clusterSize_ = cluster_size;
452 longRowThreshold_ = long_row_threshold;
453 is_matrix_structurally_symmetric_ = is_matrix_structurally_symmetric;
454 ifpack2_dump_matrix_ = ifpack2_dump_matrix;
455 localSmoothingIndices_ = localSmoothingIndices;
457 NumInnerSweeps_ = numInnerSweeps;
458 NumOuterSweeps_ = numOuterSweeps;
459 InnerSpTrsv_ = innerSpTrsv;
460 InnerDampingFactor_ = innerDampingFactor;
461 CompactForm_ = compactForm;
462 TimerForApply_ = timer_for_apply;
466template<
class MatrixType>
471 this->setParametersImpl(
const_cast<Teuchos::ParameterList&
>(pl));
475template<
class MatrixType>
476Teuchos::RCP<const Teuchos::Comm<int> >
478 TEUCHOS_TEST_FOR_EXCEPTION(
479 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::getComm: "
480 "The input matrix A is null. Please call setMatrix() with a nonnull "
481 "input matrix before calling this method.");
482 return A_->getRowMap ()->getComm ();
486template<
class MatrixType>
487Teuchos::RCP<const typename Relaxation<MatrixType>::row_matrix_type>
493template<
class MatrixType>
494Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
495 typename MatrixType::global_ordinal_type,
496 typename MatrixType::node_type> >
498 TEUCHOS_TEST_FOR_EXCEPTION(
499 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::getDomainMap: "
500 "The input matrix A is null. Please call setMatrix() with a nonnull "
501 "input matrix before calling this method.");
502 return A_->getDomainMap ();
506template<
class MatrixType>
507Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
508 typename MatrixType::global_ordinal_type,
509 typename MatrixType::node_type> >
511 TEUCHOS_TEST_FOR_EXCEPTION(
512 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::getRangeMap: "
513 "The input matrix A is null. Please call setMatrix() with a nonnull "
514 "input matrix before calling this method.");
515 return A_->getRangeMap ();
519template<
class MatrixType>
525template<
class MatrixType>
527 return(NumInitialize_);
531template<
class MatrixType>
537template<
class MatrixType>
543template<
class MatrixType>
545 return(InitializeTime_);
549template<
class MatrixType>
551 return(ComputeTime_);
555template<
class MatrixType>
561template<
class MatrixType>
563 return(ComputeFlops_);
567template<
class MatrixType>
574template<
class MatrixType>
576 TEUCHOS_TEST_FOR_EXCEPTION(
577 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::getNodeSmootherComplexity: "
578 "The input matrix A is null. Please call setMatrix() with a nonnull "
579 "input matrix, then call compute(), before calling this method.");
581 return A_->getLocalNumRows() + A_->getLocalNumEntries();
585template<
class MatrixType>
588apply (
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
589 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
597 using Teuchos::rcpFromRef;
600 TEUCHOS_TEST_FOR_EXCEPTION(
601 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::apply: "
602 "The input matrix A is null. Please call setMatrix() with a nonnull "
603 "input matrix, then call compute(), before calling this method.");
604 TEUCHOS_TEST_FOR_EXCEPTION(
607 "Ifpack2::Relaxation::apply: You must call compute() on this Ifpack2 "
608 "preconditioner instance before you may call apply(). You may call "
609 "isComputed() to find out if compute() has been called already.");
610 TEUCHOS_TEST_FOR_EXCEPTION(
611 X.getNumVectors() != Y.getNumVectors(),
613 "Ifpack2::Relaxation::apply: X and Y have different numbers of columns. "
614 "X has " << X.getNumVectors() <<
" columns, but Y has "
615 << Y.getNumVectors() <<
" columns.");
616 TEUCHOS_TEST_FOR_EXCEPTION(
617 beta != STS::zero (), std::logic_error,
618 "Ifpack2::Relaxation::apply: beta = " << beta <<
" != 0 case not "
621 Teuchos::RCP<Teuchos::Time> timer;
622 const std::string timerName (
"Ifpack2::Relaxation::apply");
623 if (TimerForApply_) {
624 timer = Teuchos::TimeMonitor::lookupCounter (timerName);
625 if (timer.is_null ()) {
626 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
630 Teuchos::Time time = Teuchos::Time(timerName);
631 double startTime = time.wallTime();
634 Teuchos::RCP<Teuchos::TimeMonitor> timeMon;
636 timeMon = Teuchos::rcp(
new Teuchos::TimeMonitor(*timer));
639 if (alpha == STS::zero ()) {
641 Y.putScalar (STS::zero ());
649 Xcopy = rcp (
new MV (X, Teuchos::Copy));
651 Xcopy = rcpFromRef (X);
659 case Ifpack2::Details::JACOBI:
660 ApplyInverseJacobi(*Xcopy,Y);
662 case Ifpack2::Details::GS:
663 ApplyInverseSerialGS(*Xcopy, Y, DoBackwardGS_ ? Tpetra::Backward : Tpetra::Forward);
665 case Ifpack2::Details::SGS:
666 ApplyInverseSerialGS(*Xcopy, Y, Tpetra::Symmetric);
668 case Ifpack2::Details::MTGS:
669 case Ifpack2::Details::GS2:
670 ApplyInverseMTGS_CrsMatrix(*Xcopy, Y, DoBackwardGS_ ? Tpetra::Backward : Tpetra::Forward);
672 case Ifpack2::Details::MTSGS:
673 case Ifpack2::Details::SGS2:
674 ApplyInverseMTGS_CrsMatrix(*Xcopy, Y, Tpetra::Symmetric);
676 case Ifpack2::Details::RICHARDSON:
677 ApplyInverseRichardson(*Xcopy,Y);
681 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
682 "Ifpack2::Relaxation::apply: Invalid preconditioner type enum value "
683 << PrecType_ <<
". Please report this bug to the Ifpack2 developers.");
685 if (alpha != STS::one ()) {
687 const double numGlobalRows = as<double> (A_->getGlobalNumRows ());
688 const double numVectors = as<double> (Y.getNumVectors ());
689 ApplyFlops_ += numGlobalRows * numVectors;
693 ApplyTime_ += (time.wallTime() - startTime);
698template<
class MatrixType>
701applyMat (
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
702 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
703 Teuchos::ETransp mode)
const
705 TEUCHOS_TEST_FOR_EXCEPTION(
706 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::applyMat: "
707 "The input matrix A is null. Please call setMatrix() with a nonnull "
708 "input matrix, then call compute(), before calling this method.");
709 TEUCHOS_TEST_FOR_EXCEPTION(
710 !
isComputed (), std::runtime_error,
"Ifpack2::Relaxation::applyMat: "
711 "isComputed() must be true before you may call applyMat(). "
712 "Please call compute() before calling this method.");
713 TEUCHOS_TEST_FOR_EXCEPTION(
714 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
715 "Ifpack2::Relaxation::applyMat: X.getNumVectors() = " << X.getNumVectors ()
716 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
717 A_->apply (X, Y, mode);
721template<
class MatrixType>
724 const char methodName[] =
"Ifpack2::Relaxation::initialize";
726 TEUCHOS_TEST_FOR_EXCEPTION
727 (A_.is_null (), std::runtime_error, methodName <<
": The "
728 "input matrix A is null. Please call setMatrix() with "
729 "a nonnull input matrix before calling this method.");
731 Teuchos::RCP<Teuchos::Time> timer =
732 Teuchos::TimeMonitor::getNewCounter (methodName);
734 double startTime = timer->wallTime();
737 Teuchos::TimeMonitor timeMon (*timer);
738 isInitialized_ =
false;
741 auto rowMap = A_->getRowMap ();
742 auto comm = rowMap.is_null () ? Teuchos::null : rowMap->getComm ();
743 IsParallel_ = ! comm.is_null () && comm->getSize () > 1;
753 Importer_ = IsParallel_ ? A_->getGraph ()->getImporter () :
757 Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
758 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A_);
759 hasBlockCrsMatrix_ = ! A_bcrs.is_null ();
762 serialGaussSeidel_ = Teuchos::null;
764 if (PrecType_ == Details::MTGS || PrecType_ == Details::MTSGS ||
765 PrecType_ == Details::GS2 || PrecType_ == Details::SGS2) {
766 auto crsMat = Details::getCrsMatrix(A_);
767 TEUCHOS_TEST_FOR_EXCEPTION
768 (crsMat.is_null(), std::logic_error, methodName <<
": "
769 "Multithreaded Gauss-Seidel methods currently only work "
770 "when the input matrix is a Tpetra::CrsMatrix.");
775 if (ifpack2_dump_matrix_) {
776 static int sequence_number = 0;
777 const std::string file_name =
"Ifpack2_MT_GS_" +
778 std::to_string (sequence_number++) +
".mtx";
779 using writer_type = Tpetra::MatrixMarket::Writer<crs_matrix_type>;
780 writer_type::writeSparseFile (file_name, crsMat);
783 this->mtKernelHandle_ = Teuchos::rcp (
new mt_kernel_handle_type ());
784 if (mtKernelHandle_->get_gs_handle () ==
nullptr) {
785 if (PrecType_ == Details::GS2 || PrecType_ == Details::SGS2)
786 mtKernelHandle_->create_gs_handle (KokkosSparse::GS_TWOSTAGE);
787 else if(this->clusterSize_ == 1) {
788 mtKernelHandle_->create_gs_handle (KokkosSparse::GS_DEFAULT, this->mtColoringAlgorithm_);
789 mtKernelHandle_->get_point_gs_handle()->set_long_row_threshold(longRowThreshold_);
792 mtKernelHandle_->create_gs_handle (KokkosSparse::CLUSTER_DEFAULT, this->clusterSize_, this->mtColoringAlgorithm_);
794 local_matrix_device_type kcsr = crsMat->getLocalMatrixDevice ();
795 if (PrecType_ == Details::GS2 || PrecType_ == Details::SGS2) {
797 mtKernelHandle_->set_gs_set_num_inner_sweeps (NumInnerSweeps_);
798 mtKernelHandle_->set_gs_set_num_outer_sweeps (NumOuterSweeps_);
799 mtKernelHandle_->set_gs_set_inner_damp_factor (InnerDampingFactor_);
800 mtKernelHandle_->set_gs_twostage (!InnerSpTrsv_, A_->getLocalNumRows ());
801 mtKernelHandle_->set_gs_twostage_compact_form (CompactForm_);
804 KokkosSparse::Experimental::gauss_seidel_symbolic(
805 mtKernelHandle_.getRawPtr (),
806 A_->getLocalNumRows (),
807 A_->getLocalNumCols (),
810 is_matrix_structurally_symmetric_);
814 InitializeTime_ += (timer->wallTime() - startTime);
816 isInitialized_ =
true;
820template <
typename BlockDiagView>
821struct InvertDiagBlocks {
822 typedef typename BlockDiagView::size_type Size;
825 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
826 template <
typename View>
827 using UnmanagedView = Kokkos::View<
typename View::data_type,
typename View::array_layout,
828 typename View::device_type, Unmanaged>;
830 typedef typename BlockDiagView::non_const_value_type Scalar;
831 typedef typename BlockDiagView::device_type Device;
832 typedef Kokkos::View<Scalar**, Kokkos::LayoutRight, Device> RWrk;
833 typedef Kokkos::View<int**, Kokkos::LayoutRight, Device> IWrk;
835 UnmanagedView<BlockDiagView> block_diag_;
839 UnmanagedView<RWrk> rwrk_;
841 UnmanagedView<IWrk> iwrk_;
844 InvertDiagBlocks (BlockDiagView& block_diag)
845 : block_diag_(block_diag)
847 const auto blksz = block_diag.extent(1);
848 Kokkos::resize(rwrk_buf_, block_diag_.extent(0), blksz);
850 Kokkos::resize(iwrk_buf_, block_diag_.extent(0), blksz);
854 KOKKOS_INLINE_FUNCTION
855 void operator() (
const Size i,
int& jinfo)
const {
856 auto D_cur = Kokkos::subview(block_diag_, i, Kokkos::ALL(), Kokkos::ALL());
857 auto ipiv = Kokkos::subview(iwrk_, i, Kokkos::ALL());
858 auto work = Kokkos::subview(rwrk_, i, Kokkos::ALL());
860 Tpetra::GETF2(D_cur, ipiv, info);
865 Tpetra::GETRI(D_cur, ipiv, work, info);
871template<
class MatrixType>
872void Relaxation<MatrixType>::computeBlockCrs ()
875 using Teuchos::Array;
876 using Teuchos::ArrayRCP;
877 using Teuchos::ArrayView;
882 using Teuchos::REDUCE_MAX;
883 using Teuchos::REDUCE_MIN;
884 using Teuchos::REDUCE_SUM;
885 using Teuchos::rcp_dynamic_cast;
886 using Teuchos::reduceAll;
889 const std::string timerName (
"Ifpack2::Relaxation::computeBlockCrs");
890 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
891 if (timer.is_null ()) {
892 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
894 double startTime = timer->wallTime();
896 Teuchos::TimeMonitor timeMon (*timer);
898 TEUCHOS_TEST_FOR_EXCEPTION(
899 A_.is_null (), std::runtime_error,
"Ifpack2::Relaxation::"
900 "computeBlockCrs: The input matrix A is null. Please call setMatrix() "
901 "with a nonnull input matrix, then call initialize(), before calling "
903 auto blockCrsA = rcp_dynamic_cast<const block_crs_matrix_type> (A_);
904 TEUCHOS_TEST_FOR_EXCEPTION(
905 blockCrsA.is_null(), std::logic_error,
"Ifpack2::Relaxation::"
906 "computeBlockCrs: A_ is not a BlockCrsMatrix, but it should be if we "
907 "got this far. Please report this bug to the Ifpack2 developers.");
909 const scalar_type one = STS::one ();
914 const LO lclNumMeshRows =
915 blockCrsA->getCrsGraph ().getLocalNumRows ();
916 const LO blockSize = blockCrsA->getBlockSize ();
918 if (! savedDiagOffsets_) {
919 blockDiag_ = block_diag_type ();
920 blockDiag_ = block_diag_type (
"Ifpack2::Relaxation::blockDiag_",
921 lclNumMeshRows, blockSize, blockSize);
922 if (Teuchos::as<LO>(diagOffsets_.extent (0) ) < lclNumMeshRows) {
925 diagOffsets_ = Kokkos::View<size_t*, device_type> ();
926 diagOffsets_ = Kokkos::View<size_t*, device_type> (
"offsets", lclNumMeshRows);
928 blockCrsA->getCrsGraph ().getLocalDiagOffsets (diagOffsets_);
929 TEUCHOS_TEST_FOR_EXCEPTION
930 (
static_cast<size_t> (diagOffsets_.extent (0)) !=
931 static_cast<size_t> (blockDiag_.extent (0)),
932 std::logic_error,
"diagOffsets_.extent(0) = " <<
933 diagOffsets_.extent (0) <<
" != blockDiag_.extent(0) = "
934 << blockDiag_.extent (0) <<
935 ". Please report this bug to the Ifpack2 developers.");
936 savedDiagOffsets_ =
true;
938 blockCrsA->getLocalDiagCopy (blockDiag_, diagOffsets_);
945 unmanaged_block_diag_type blockDiag = blockDiag_;
947 if (DoL1Method_ && IsParallel_) {
948 const scalar_type two = one + one;
949 const size_t maxLength = A_->getLocalMaxNumRowEntries ();
950 nonconst_local_inds_host_view_type indices (
"indices",maxLength);
951 nonconst_values_host_view_type values_ (
"values",maxLength * blockSize * blockSize);
952 size_t numEntries = 0;
954 for (LO i = 0; i < lclNumMeshRows; ++i) {
956 blockCrsA->getLocalRowCopy (i, indices, values_, numEntries);
957 scalar_type * values =
reinterpret_cast<scalar_type*
>(values_.data());
959 auto diagBlock = Kokkos::subview (blockDiag, i, ALL (), ALL ());
960 for (LO subRow = 0; subRow < blockSize; ++subRow) {
961 magnitude_type diagonal_boost = STM::zero ();
962 for (
size_t k = 0 ; k < numEntries ; ++k) {
963 if (indices[k] >= lclNumMeshRows) {
964 const size_t offset = blockSize*blockSize*k + subRow*blockSize;
965 for (LO subCol = 0; subCol < blockSize; ++subCol) {
966 diagonal_boost += STS::magnitude (values[offset+subCol] / two);
970 if (STS::magnitude (diagBlock(subRow, subRow)) < L1Eta_ * diagonal_boost) {
971 diagBlock(subRow, subRow) += diagonal_boost;
979 Impl::InvertDiagBlocks<unmanaged_block_diag_type> idb(blockDiag);
980 typedef typename unmanaged_block_diag_type::execution_space exec_space;
981 typedef Kokkos::RangePolicy<exec_space, LO> range_type;
983 Kokkos::parallel_reduce (range_type (0, lclNumMeshRows), idb, info);
988 TEUCHOS_TEST_FOR_EXCEPTION
989 (info > 0, std::runtime_error,
"GETF2 or GETRI failed on = " << info
990 <<
" diagonal blocks.");
995#ifdef HAVE_IFPACK2_DEBUG
996 const int numResults = 2;
998 int lclResults[2], gblResults[2];
999 lclResults[0] = info;
1000 lclResults[1] = -info;
1003 reduceAll<int, int> (* (A_->getGraph ()->getComm ()), REDUCE_MIN,
1004 numResults, lclResults, gblResults);
1005 TEUCHOS_TEST_FOR_EXCEPTION
1006 (gblResults[0] != 0 || gblResults[1] != 0, std::runtime_error,
1007 "Ifpack2::Relaxation::compute: When processing the input "
1008 "Tpetra::BlockCrsMatrix, one or more diagonal block LU factorizations "
1009 "failed on one or more (MPI) processes.");
1011 serialGaussSeidel_ = rcp(
new SerialGaussSeidel(blockCrsA, blockDiag_, localSmoothingIndices_, DampingFactor_));
1014 ComputeTime_ += (timer->wallTime() - startTime);
1019template<
class MatrixType>
1022 using Teuchos::Array;
1023 using Teuchos::ArrayRCP;
1024 using Teuchos::ArrayView;
1026 using Teuchos::Comm;
1029 using Teuchos::REDUCE_MAX;
1030 using Teuchos::REDUCE_MIN;
1031 using Teuchos::REDUCE_SUM;
1032 using Teuchos::reduceAll;
1036 using IST =
typename vector_type::impl_scalar_type;
1037 using KAT = Kokkos::ArithTraits<IST>;
1039 const char methodName[] =
"Ifpack2::Relaxation::compute";
1043 TEUCHOS_TEST_FOR_EXCEPTION
1044 (A_.is_null (), std::runtime_error, methodName <<
": "
1045 "The input matrix A is null. Please call setMatrix() with a nonnull "
1046 "input matrix, then call initialize(), before calling this method.");
1053 if (hasBlockCrsMatrix_) {
1058 Teuchos::RCP<Teuchos::Time> timer =
1059 Teuchos::TimeMonitor::getNewCounter (methodName);
1060 double startTime = timer->wallTime();
1063 Teuchos::TimeMonitor timeMon (*timer);
1072 IST oneOverMinDiagVal = KAT::one () /
static_cast<IST
> (SmallTraits<scalar_type>::eps ());
1073 if ( MinDiagonalValue_ != zero)
1074 oneOverMinDiagVal = KAT::one () /
static_cast<IST
> (MinDiagonalValue_);
1077 const magnitude_type minDiagValMag = STS::magnitude (MinDiagonalValue_);
1079 const LO numMyRows =
static_cast<LO
> (A_->getLocalNumRows ());
1081 TEUCHOS_TEST_FOR_EXCEPTION
1082 (NumSweeps_ < 0, std::logic_error, methodName
1083 <<
": NumSweeps_ = " << NumSweeps_ <<
" < 0. "
1084 "Please report this bug to the Ifpack2 developers.");
1085 IsComputed_ =
false;
1087 if (Diagonal_.is_null()) {
1090 Diagonal_ = rcp (
new vector_type (A_->getRowMap (),
false));
1093 if (checkDiagEntries_) {
1099 size_t numSmallDiagEntries = 0;
1100 size_t numZeroDiagEntries = 0;
1101 size_t numNegDiagEntries = 0;
1108 A_->getLocalDiagCopy (*Diagonal_);
1109 std::unique_ptr<vector_type> origDiag;
1110 origDiag = std::unique_ptr<vector_type> (
new vector_type (*Diagonal_, Teuchos::Copy));
1111 auto diag2d = Diagonal_->getLocalViewHost(Tpetra::Access::ReadWrite);
1112 auto diag = Kokkos::subview(diag2d, Kokkos::ALL(), 0);
1119 if (numMyRows != 0) {
1121 minMagDiagEntryMag = d_0_mag;
1122 maxMagDiagEntryMag = d_0_mag;
1131 for (LO i = 0; i < numMyRows; ++i) {
1132 const IST d_i = diag(i);
1136 const auto d_i_real = getRealValue (d_i);
1140 if (d_i_real < STM::zero ()) {
1141 ++numNegDiagEntries;
1143 if (d_i_mag < minMagDiagEntryMag) {
1144 minMagDiagEntryMag = d_i_mag;
1146 if (d_i_mag > maxMagDiagEntryMag) {
1147 maxMagDiagEntryMag = d_i_mag;
1150 if (fixTinyDiagEntries_) {
1152 if (d_i_mag <= minDiagValMag) {
1153 ++numSmallDiagEntries;
1154 if (d_i_mag == STM::zero ()) {
1155 ++numZeroDiagEntries;
1157 diag(i) = oneOverMinDiagVal;
1160 diag(i) = KAT::one () / d_i;
1165 if (d_i_mag <= minDiagValMag) {
1166 ++numSmallDiagEntries;
1167 if (d_i_mag == STM::zero ()) {
1168 ++numZeroDiagEntries;
1171 diag(i) = KAT::one () / d_i;
1189 const Comm<int>& comm = * (A_->getRowMap ()->
getComm ());
1192 Array<magnitude_type> localVals (2);
1193 localVals[0] = minMagDiagEntryMag;
1195 localVals[1] = -maxMagDiagEntryMag;
1196 Array<magnitude_type> globalVals (2);
1197 reduceAll<int, magnitude_type> (comm, REDUCE_MIN, 2,
1198 localVals.getRawPtr (),
1199 globalVals.getRawPtr ());
1200 globalMinMagDiagEntryMag_ = globalVals[0];
1201 globalMaxMagDiagEntryMag_ = -globalVals[1];
1203 Array<size_t> localCounts (3);
1204 localCounts[0] = numSmallDiagEntries;
1205 localCounts[1] = numZeroDiagEntries;
1206 localCounts[2] = numNegDiagEntries;
1207 Array<size_t> globalCounts (3);
1208 reduceAll<int, size_t> (comm, REDUCE_SUM, 3,
1209 localCounts.getRawPtr (),
1210 globalCounts.getRawPtr ());
1211 globalNumSmallDiagEntries_ = globalCounts[0];
1212 globalNumZeroDiagEntries_ = globalCounts[1];
1213 globalNumNegDiagEntries_ = globalCounts[2];
1218 diff.reciprocal (*origDiag);
1219 diff.update (-one, *Diagonal_, one);
1220 globalDiagNormDiff_ = diff.norm2 ();
1227 bool debugAgainstSlowPath =
false;
1229 auto crsMat = Details::getCrsMatrix(A_);
1231 if (crsMat.get() && crsMat->isFillComplete ()) {
1236 if (invDiagKernel_.is_null())
1239 invDiagKernel_->setMatrix(crsMat);
1240 invDiagKernel_->compute(*Diagonal_,
1241 DoL1Method_ && IsParallel_, L1Eta_,
1242 fixTinyDiagEntries_, minDiagValMag);
1243 savedDiagOffsets_ =
true;
1247#ifdef HAVE_IFPACK2_DEBUG
1248 debugAgainstSlowPath =
true;
1252 if (crsMat.is_null() || ! crsMat->isFillComplete () || debugAgainstSlowPath) {
1260 if (debugAgainstSlowPath)
1275 if (DoL1Method_ && IsParallel_) {
1277 auto diag =
Diagonal->getLocalViewHost(Tpetra::Access::ReadWrite);
1279 const size_t maxLength = A_row.getLocalMaxNumRowEntries ();
1280 nonconst_local_inds_host_view_type indices(
"indices",maxLength);
1281 nonconst_values_host_view_type values(
"values",maxLength);
1284 for (LO i = 0; i < numMyRows; ++i) {
1285 A_row.getLocalRowCopy (i, indices, values, numEntries);
1287 for (
size_t k = 0 ; k < numEntries; ++k) {
1288 if (indices[k] >= numMyRows) {
1289 diagonal_boost += STS::magnitude (values[k] / two);
1292 if (KAT::magnitude (diag(i, 0)) < L1Eta_ * diagonal_boost) {
1293 diag(i, 0) += diagonal_boost;
1301 if (fixTinyDiagEntries_) {
1305 auto localDiag =
Diagonal->getLocalViewDevice(Tpetra::Access::ReadWrite);
1306 Kokkos::parallel_for(Kokkos::RangePolicy<MyExecSpace>(0, localDiag.extent(0)),
1308 auto d_i = localDiag(i, 0);
1311 if (d_i_mag <= minDiagValMag) {
1312 d_i = oneOverMinDiagVal;
1317 d_i = IST (KAT::one () / d_i);
1319 localDiag(i, 0) = d_i;
1326 if (debugAgainstSlowPath) {
1328 Diagonal->update (STS::one (), *Diagonal_, -STS::one ());
1332 TEUCHOS_TEST_FOR_EXCEPTION
1333 (err > 100*STM::eps(), std::logic_error, methodName <<
": "
1334 <<
"\"fast-path\" diagonal computation failed. "
1335 "\\|D1 - D2\\|_inf = " << err <<
".");
1342 if (STS::isComplex) {
1343 ComputeFlops_ += 4.0 * numMyRows;
1346 ComputeFlops_ += numMyRows;
1349 if (PrecType_ == Ifpack2::Details::MTGS ||
1350 PrecType_ == Ifpack2::Details::MTSGS ||
1351 PrecType_ == Ifpack2::Details::GS2 ||
1352 PrecType_ == Ifpack2::Details::SGS2) {
1355 using scalar_view_t =
typename local_matrix_device_type::values_type;
1357 TEUCHOS_TEST_FOR_EXCEPTION
1358 (crsMat.is_null(), std::logic_error, methodName <<
": "
1359 "Multithreaded Gauss-Seidel methods currently only work "
1360 "when the input matrix is a Tpetra::CrsMatrix.");
1361 local_matrix_device_type kcsr = crsMat->getLocalMatrixDevice ();
1365 auto diagView_2d = Diagonal_->getLocalViewDevice (Tpetra::Access::ReadWrite);
1366 scalar_view_t diagView_1d = Kokkos::subview (diagView_2d, Kokkos::ALL (), 0);
1367 KokkosSparse::Experimental::gauss_seidel_numeric(
1368 mtKernelHandle_.getRawPtr (),
1369 A_->getLocalNumRows (),
1370 A_->getLocalNumCols (),
1375 is_matrix_structurally_symmetric_);
1377 else if(PrecType_ == Ifpack2::Details::GS || PrecType_ == Ifpack2::Details::SGS) {
1379 serialGaussSeidel_ = rcp(
new SerialGaussSeidel(crsMat, Diagonal_, localSmoothingIndices_, DampingFactor_));
1381 serialGaussSeidel_ = rcp(
new SerialGaussSeidel(A_, Diagonal_, localSmoothingIndices_, DampingFactor_));
1385 ComputeTime_ += (timer->wallTime() - startTime);
1391template<
class MatrixType>
1394ApplyInverseRichardson (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1395 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const
1398 const double numGlobalRows = as<double> (A_->getGlobalNumRows ());
1399 const double numVectors = as<double> (X.getNumVectors ());
1400 if (ZeroStartingSolution_) {
1404 Y.scale(DampingFactor_,X);
1410 double flopUpdate = 0.0;
1411 if (DampingFactor_ == STS::one ()) {
1413 flopUpdate = numGlobalRows * numVectors;
1417 flopUpdate = numGlobalRows * numVectors;
1419 ApplyFlops_ += flopUpdate;
1420 if (NumSweeps_ == 1) {
1426 const int startSweep = ZeroStartingSolution_ ? 1 : 0;
1429 updateCachedMultiVector(Y.getMap(),as<size_t>(numVectors));
1431 for (
int j = startSweep; j < NumSweeps_; ++j) {
1433 Tpetra::Details::residual(*A_,Y,X,*cachedMV_);
1434 Y.update(DampingFactor_,*cachedMV_,STS::one());
1448 const double numGlobalNonzeros = as<double> (A_->getGlobalNumEntries ());
1449 const double dampingFlops = (DampingFactor_ == STS::one ()) ? 0.0 : 1.0;
1450 ApplyFlops_ += as<double> (NumSweeps_ - startSweep) * numVectors *
1451 (2.0 * numGlobalNonzeros + dampingFlops);
1455template<
class MatrixType>
1458ApplyInverseJacobi (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1459 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const
1462 if (hasBlockCrsMatrix_) {
1463 ApplyInverseJacobi_BlockCrsMatrix (X, Y);
1467 const double numGlobalRows = as<double> (A_->getGlobalNumRows ());
1468 const double numVectors = as<double> (X.getNumVectors ());
1469 if (ZeroStartingSolution_) {
1474 Y.elementWiseMultiply (DampingFactor_, *Diagonal_, X, STS::zero ());
1480 double flopUpdate = 0.0;
1481 if (DampingFactor_ == STS::one ()) {
1483 flopUpdate = numGlobalRows * numVectors;
1487 flopUpdate = 2.0 * numGlobalRows * numVectors;
1489 ApplyFlops_ += flopUpdate;
1490 if (NumSweeps_ == 1) {
1496 const int startSweep = ZeroStartingSolution_ ? 1 : 0;
1499 updateCachedMultiVector(Y.getMap(),as<size_t>(numVectors));
1501 for (
int j = startSweep; j < NumSweeps_; ++j) {
1503 Tpetra::Details::residual(*A_,Y,X,*cachedMV_);
1504 Y.elementWiseMultiply (DampingFactor_, *Diagonal_, *cachedMV_, STS::one ());
1518 const double numGlobalNonzeros = as<double> (A_->getGlobalNumEntries ());
1519 const double dampingFlops = (DampingFactor_ == STS::one ()) ? 0.0 : 1.0;
1520 ApplyFlops_ += as<double> (NumSweeps_ - startSweep) * numVectors *
1521 (2.0 * numGlobalRows + 2.0 * numGlobalNonzeros + dampingFlops);
1524template<
class MatrixType>
1531 Tpetra::MultiVector<scalar_type,
1536 using Tpetra::BlockMultiVector;
1540 const block_crs_matrix_type* blockMatConst =
1541 dynamic_cast<const block_crs_matrix_type*
> (A_.getRawPtr ());
1542 TEUCHOS_TEST_FOR_EXCEPTION
1543 (blockMatConst ==
nullptr, std::logic_error,
"This method should "
1544 "never be called if the matrix A_ is not a BlockCrsMatrix. "
1545 "Please report this bug to the Ifpack2 developers.");
1550 block_crs_matrix_type* blockMat =
1551 const_cast<block_crs_matrix_type*
> (blockMatConst);
1553 auto meshRowMap = blockMat->getRowMap ();
1554 auto meshColMap = blockMat->getColMap ();
1557 BMV X_blk (X, *meshColMap, blockSize);
1558 BMV Y_blk (Y, *meshRowMap, blockSize);
1560 if (ZeroStartingSolution_) {
1564 Y_blk.blockWiseMultiply (DampingFactor_, blockDiag_, X_blk);
1565 if (NumSweeps_ == 1) {
1570 auto pointRowMap = Y.getMap ();
1571 const size_t numVecs = X.getNumVectors ();
1575 BMV A_times_Y (*meshRowMap, *pointRowMap, blockSize, numVecs);
1579 const int startSweep = ZeroStartingSolution_ ? 1 : 0;
1581 for (
int j = startSweep; j < NumSweeps_; ++j) {
1582 blockMat->applyBlock (Y_blk, A_times_Y);
1584 Y_blk.blockJacobiUpdate (DampingFactor_, blockDiag_,
1585 X_blk, A_times_Y, STS::one ());
1589template<
class MatrixType>
1592ApplyInverseSerialGS (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1593 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1594 Tpetra::ESweepDirection direction)
const
1606 auto blockCrsMat = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A_);
1607 auto crsMat = Details::getCrsMatrix(A_);
1608 if (blockCrsMat.get()) {
1609 const_cast<this_type&
> (*this).ApplyInverseSerialGS_BlockCrsMatrix (*blockCrsMat, X, Y, direction);
1611 else if (crsMat.get()) {
1612 ApplyInverseSerialGS_CrsMatrix (*crsMat, X, Y, direction);
1615 ApplyInverseSerialGS_RowMatrix (X, Y, direction);
1620template<
class MatrixType>
1624 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1625 Tpetra::ESweepDirection direction)
const {
1626 using Teuchos::Array;
1627 using Teuchos::ArrayRCP;
1628 using Teuchos::ArrayView;
1632 using Teuchos::rcpFromRef;
1633 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
1638 if (ZeroStartingSolution_) {
1639 Y.putScalar (STS::zero ());
1642 size_t NumVectors = X.getNumVectors();
1646 if (Importer_.is_null ()) {
1647 updateCachedMultiVector (Y.getMap (), NumVectors);
1650 updateCachedMultiVector (Importer_->getTargetMap (), NumVectors);
1655 Y2 = rcpFromRef (Y);
1658 for (
int j = 0; j < NumSweeps_; ++j) {
1661 if (Importer_.is_null ()) {
1667 Y2->doImport (Y, *Importer_, Tpetra::INSERT);
1670 serialGaussSeidel_->apply(*Y2, X, direction);
1674 Tpetra::deep_copy (Y, *Y2);
1679 const double dampingFlops = (DampingFactor_ == STS::one()) ? 0.0 : 1.0;
1680 const double numVectors = as<double> (X.getNumVectors ());
1681 const double numGlobalRows = as<double> (A_->getGlobalNumRows ());
1682 const double numGlobalNonzeros = as<double> (A_->getGlobalNumEntries ());
1683 ApplyFlops_ += 2.0 * NumSweeps_ * numVectors *
1684 (2.0 * numGlobalRows + 2.0 * numGlobalNonzeros + dampingFlops);
1688template<
class MatrixType>
1692 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
1693 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1694 Tpetra::ESweepDirection direction)
const
1696 using Teuchos::null;
1699 using Teuchos::rcpFromRef;
1700 using Teuchos::rcp_const_cast;
1701 typedef scalar_type Scalar;
1702 const char prefix[] =
"Ifpack2::Relaxation::SerialGS: ";
1703 const scalar_type ZERO = Teuchos::ScalarTraits<Scalar>::zero ();
1705 TEUCHOS_TEST_FOR_EXCEPTION(
1706 ! A.isFillComplete (), std::runtime_error,
1707 prefix <<
"The matrix is not fill complete.");
1708 TEUCHOS_TEST_FOR_EXCEPTION(
1709 NumSweeps_ < 0, std::invalid_argument,
1710 prefix <<
"The number of sweeps must be nonnegative, "
1711 "but you provided numSweeps = " << NumSweeps_ <<
" < 0.");
1713 if (NumSweeps_ == 0) {
1717 RCP<const import_type> importer = A.getGraph ()->getImporter ();
1719 RCP<const map_type> domainMap = A.getDomainMap ();
1720 RCP<const map_type> rangeMap = A.getRangeMap ();
1721 RCP<const map_type> rowMap = A.getGraph ()->getRowMap ();
1722 RCP<const map_type> colMap = A.getGraph ()->getColMap ();
1724#ifdef HAVE_IFPACK2_DEBUG
1729 TEUCHOS_TEST_FOR_EXCEPTION(
1730 ! X.getMap ()->isSameAs (*domainMap), std::runtime_error,
1731 "Tpetra::CrsMatrix::gaussSeidelCopy requires that the input "
1732 "multivector X be in the domain Map of the matrix.");
1733 TEUCHOS_TEST_FOR_EXCEPTION(
1734 ! B.getMap ()->isSameAs (*rangeMap), std::runtime_error,
1735 "Tpetra::CrsMatrix::gaussSeidelCopy requires that the input "
1736 "B be in the range Map of the matrix.");
1737 TEUCHOS_TEST_FOR_EXCEPTION(
1738 ! Diagonal_->getMap ()->isSameAs (*rowMap), std::runtime_error,
1739 "Tpetra::CrsMatrix::gaussSeidelCopy requires that the input "
1740 "D be in the row Map of the matrix.");
1741 TEUCHOS_TEST_FOR_EXCEPTION(
1742 ! rowMap->isSameAs (*rangeMap), std::runtime_error,
1743 "Tpetra::CrsMatrix::gaussSeidelCopy requires that the row Map and the "
1744 "range Map be the same (in the sense of Tpetra::Map::isSameAs).");
1745 TEUCHOS_TEST_FOR_EXCEPTION(
1746 ! domainMap->isSameAs (*rangeMap), std::runtime_error,
1747 "Tpetra::CrsMatrix::gaussSeidelCopy requires that the domain Map and "
1748 "the range Map of the matrix be the same.");
1760 RCP<multivector_type> X_colMap;
1761 RCP<multivector_type> X_domainMap;
1762 bool copyBackOutput =
false;
1763 if (importer.is_null ()) {
1764 X_colMap = Teuchos::rcpFromRef (X);
1765 X_domainMap = Teuchos::rcpFromRef (X);
1771 if (ZeroStartingSolution_) {
1772 X_colMap->putScalar (ZERO);
1777 updateCachedMultiVector(colMap, X.getNumVectors());
1778 X_colMap = cachedMV_;
1779 X_domainMap = X_colMap->offsetViewNonConst (domainMap, 0);
1781 if (ZeroStartingSolution_) {
1783 X_colMap->putScalar (ZERO);
1793 X_colMap->doImport (X, *importer, Tpetra::INSERT);
1795 copyBackOutput =
true;
1798 for (
int sweep = 0; sweep < NumSweeps_; ++sweep) {
1799 if (! importer.is_null () && sweep > 0) {
1802 X_colMap->doImport (*X_domainMap, *importer, Tpetra::INSERT);
1805 serialGaussSeidel_->apply(*X_colMap, B, direction);
1808 if (copyBackOutput) {
1810 deep_copy (X , *X_domainMap);
1811 }
catch (std::exception& e) {
1812 TEUCHOS_TEST_FOR_EXCEPTION(
1813 true, std::runtime_error, prefix <<
"deep_copy(X, *X_domainMap) "
1814 "threw an exception: " << e.what ());
1829 const double dampingFlops = (DampingFactor_ == STS::one()) ? 0.0 : 1.0;
1830 const double numVectors = X.getNumVectors ();
1831 const double numGlobalRows = A_->getGlobalNumRows ();
1832 const double numGlobalNonzeros = A_->getGlobalNumEntries ();
1833 ApplyFlops_ += NumSweeps_ * numVectors *
1834 (2.0 * numGlobalRows + 2.0 * numGlobalNonzeros + dampingFlops);
1837template<
class MatrixType>
1841 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1842 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1843 Tpetra::ESweepDirection direction)
1845 using Tpetra::INSERT;
1848 using Teuchos::rcpFromRef;
1857 block_multivector_type yBlock(Y, *(A.getGraph ()->getDomainMap()), A.getBlockSize());
1858 const block_multivector_type xBlock(X, *(A.getColMap ()), A.getBlockSize());
1860 bool performImport =
false;
1861 RCP<block_multivector_type> yBlockCol;
1862 if (Importer_.is_null()) {
1863 yBlockCol = rcpFromRef(yBlock);
1866 if (yBlockColumnPointMap_.is_null() ||
1867 yBlockColumnPointMap_->getNumVectors() != yBlock.getNumVectors() ||
1868 yBlockColumnPointMap_->getBlockSize() != yBlock.getBlockSize()) {
1869 yBlockColumnPointMap_ =
1870 rcp(
new block_multivector_type(*(A.getColMap()), A.getBlockSize(),
1873 yBlockCol = yBlockColumnPointMap_;
1874 if (pointImporter_.is_null()) {
1875 auto srcMap = rcp(
new map_type(yBlock.getPointMap()));
1876 auto tgtMap = rcp(
new map_type(yBlockCol->getPointMap()));
1877 pointImporter_ = rcp(
new import_type(srcMap, tgtMap));
1879 performImport =
true;
1882 multivector_type yBlock_mv;
1883 multivector_type yBlockCol_mv;
1884 RCP<const multivector_type> yBlockColPointDomain;
1885 if (performImport) {
1886 yBlock_mv = yBlock.getMultiVectorView();
1887 yBlockCol_mv = yBlockCol->getMultiVectorView();
1888 yBlockColPointDomain =
1889 yBlockCol_mv.offsetView(A.getDomainMap(), 0);
1892 if (ZeroStartingSolution_) {
1893 yBlockCol->putScalar(STS::zero ());
1895 else if (performImport) {
1896 yBlockCol_mv.doImport(yBlock_mv, *pointImporter_, Tpetra::INSERT);
1899 for (
int sweep = 0; sweep < NumSweeps_; ++sweep) {
1900 if (performImport && sweep > 0) {
1901 yBlockCol_mv.doImport(yBlock_mv, *pointImporter_, Tpetra::INSERT);
1903 serialGaussSeidel_->applyBlock(*yBlockCol, xBlock, direction);
1904 if (performImport) {
1905 Tpetra::deep_copy(Y, *yBlockColPointDomain);
1910template<
class MatrixType>
1914 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
1915 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1916 Tpetra::ESweepDirection direction)
const
1918 using Teuchos::null;
1921 using Teuchos::rcpFromRef;
1922 using Teuchos::rcp_const_cast;
1925 typedef scalar_type Scalar;
1927 const char prefix[] =
"Ifpack2::Relaxation::(reordered)MTGaussSeidel: ";
1928 const Scalar ZERO = Teuchos::ScalarTraits<Scalar>::zero ();
1930 auto crsMat = Details::getCrsMatrix(A_);
1931 TEUCHOS_TEST_FOR_EXCEPTION
1932 (crsMat.is_null(), std::logic_error,
"Ifpack2::Relaxation::apply: "
1933 "Multithreaded Gauss-Seidel methods currently only work when the "
1934 "input matrix is a Tpetra::CrsMatrix.");
1937 TEUCHOS_TEST_FOR_EXCEPTION
1938 (! localSmoothingIndices_.is_null (), std::logic_error,
1939 "Ifpack2's implementation of Multithreaded Gauss-Seidel does not "
1940 "implement the case where the user supplies an iteration order. "
1941 "This error used to appear as \"MT GaussSeidel ignores the given "
1943 "I tried to add more explanation, but I didn't implement \"MT "
1944 "GaussSeidel\" [sic]. "
1945 "You'll have to ask the person who did.");
1947 TEUCHOS_TEST_FOR_EXCEPTION
1948 (! crsMat->isFillComplete (), std::runtime_error, prefix <<
"The "
1949 "input CrsMatrix is not fill complete. Please call fillComplete "
1950 "on the matrix, then do setup again, before calling apply(). "
1951 "\"Do setup\" means that if only the matrix's values have changed "
1952 "since last setup, you need only call compute(). If the matrix's "
1953 "structure may also have changed, you must first call initialize(), "
1954 "then call compute(). If you have not set up this preconditioner "
1955 "for this matrix before, you must first call initialize(), then "
1957 TEUCHOS_TEST_FOR_EXCEPTION
1958 (NumSweeps_ < 0, std::logic_error, prefix <<
": NumSweeps_ = "
1959 << NumSweeps_ <<
" < 0. Please report this bug to the Ifpack2 "
1961 if (NumSweeps_ == 0) {
1965 RCP<const import_type> importer = crsMat->getGraph ()->getImporter ();
1966 TEUCHOS_TEST_FOR_EXCEPTION(
1967 ! crsMat->getGraph ()->getExporter ().is_null (), std::runtime_error,
1968 "This method's implementation currently requires that the matrix's row, "
1969 "domain, and range Maps be the same. This cannot be the case, because "
1970 "the matrix has a nontrivial Export object.");
1972 RCP<const map_type> domainMap = crsMat->getDomainMap ();
1973 RCP<const map_type> rangeMap = crsMat->getRangeMap ();
1974 RCP<const map_type> rowMap = crsMat->getGraph ()->getRowMap ();
1975 RCP<const map_type> colMap = crsMat->getGraph ()->getColMap ();
1977#ifdef HAVE_IFPACK2_DEBUG
1982 TEUCHOS_TEST_FOR_EXCEPTION(
1983 ! X.getMap ()->isSameAs (*domainMap), std::runtime_error,
1984 "Ifpack2::Relaxation::MTGaussSeidel requires that the input "
1985 "multivector X be in the domain Map of the matrix.");
1986 TEUCHOS_TEST_FOR_EXCEPTION(
1987 ! B.getMap ()->isSameAs (*rangeMap), std::runtime_error,
1988 "Ifpack2::Relaxation::MTGaussSeidel requires that the input "
1989 "B be in the range Map of the matrix.");
1990 TEUCHOS_TEST_FOR_EXCEPTION(
1991 ! rowMap->isSameAs (*rangeMap), std::runtime_error,
1992 "Ifpack2::Relaxation::MTGaussSeidel requires that the row Map and the "
1993 "range Map be the same (in the sense of Tpetra::Map::isSameAs).");
1994 TEUCHOS_TEST_FOR_EXCEPTION(
1995 ! domainMap->isSameAs (*rangeMap), std::runtime_error,
1996 "Ifpack2::Relaxation::MTGaussSeidel requires that the domain Map and "
1997 "the range Map of the matrix be the same.");
2013 RCP<multivector_type> X_colMap;
2014 RCP<multivector_type> X_domainMap;
2015 bool copyBackOutput =
false;
2016 if (importer.is_null ()) {
2017 if (X.isConstantStride ()) {
2018 X_colMap = rcpFromRef (X);
2019 X_domainMap = rcpFromRef (X);
2026 if (ZeroStartingSolution_) {
2027 X_colMap->putScalar (ZERO);
2036 updateCachedMultiVector(colMap,as<size_t>(X.getNumVectors()));
2037 X_colMap = cachedMV_;
2041 X_domainMap = X_colMap;
2042 if (! ZeroStartingSolution_) {
2044 deep_copy (*X_domainMap , X);
2045 }
catch (std::exception& e) {
2046 std::ostringstream os;
2047 os <<
"Ifpack2::Relaxation::MTGaussSeidel: "
2048 "deep_copy(*X_domainMap, X) threw an exception: "
2049 << e.what () <<
".";
2050 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error, e.what ());
2053 copyBackOutput =
true;
2067 updateCachedMultiVector(colMap,as<size_t>(X.getNumVectors()));
2068 X_colMap = cachedMV_;
2070 X_domainMap = X_colMap->offsetViewNonConst (domainMap, 0);
2072 if (ZeroStartingSolution_) {
2074 X_colMap->putScalar (ZERO);
2084 X_colMap->doImport (X, *importer, Tpetra::CombineMode::INSERT);
2086 copyBackOutput =
true;
2092 RCP<const multivector_type> B_in;
2093 if (B.isConstantStride ()) {
2094 B_in = rcpFromRef (B);
2100 RCP<multivector_type> B_in_nonconst = rcp (
new multivector_type(rowMap, B.getNumVectors()));
2102 deep_copy (*B_in_nonconst, B);
2103 }
catch (std::exception& e) {
2104 std::ostringstream os;
2105 os <<
"Ifpack2::Relaxation::MTGaussSeidel: "
2106 "deep_copy(*B_in_nonconst, B) threw an exception: "
2107 << e.what () <<
".";
2108 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error, e.what ());
2110 B_in = rcp_const_cast<const multivector_type> (B_in_nonconst);
2124 local_matrix_device_type kcsr = crsMat->getLocalMatrixDevice ();
2126 bool update_y_vector =
true;
2128 bool zero_x_vector =
false;
2130 for (
int sweep = 0; sweep < NumSweeps_; ++sweep) {
2131 if (! importer.is_null () && sweep > 0) {
2134 X_colMap->doImport (*X_domainMap, *importer, Tpetra::CombineMode::INSERT);
2137 if (direction == Tpetra::Symmetric) {
2138 KokkosSparse::Experimental::symmetric_gauss_seidel_apply
2139 (mtKernelHandle_.getRawPtr(), A_->getLocalNumRows(), A_->getLocalNumCols(),
2140 kcsr.graph.row_map, kcsr.graph.entries, kcsr.values,
2141 X_colMap->getLocalViewDevice(Tpetra::Access::ReadWrite),
2142 B_in->getLocalViewDevice(Tpetra::Access::ReadOnly),
2143 zero_x_vector, update_y_vector, DampingFactor_, 1);
2145 else if (direction == Tpetra::Forward) {
2146 KokkosSparse::Experimental::forward_sweep_gauss_seidel_apply
2147 (mtKernelHandle_.getRawPtr(), A_->getLocalNumRows(), A_->getLocalNumCols(),
2148 kcsr.graph.row_map,kcsr.graph.entries, kcsr.values,
2149 X_colMap->getLocalViewDevice(Tpetra::Access::ReadWrite),
2150 B_in->getLocalViewDevice(Tpetra::Access::ReadOnly),
2151 zero_x_vector, update_y_vector, DampingFactor_, 1);
2153 else if (direction == Tpetra::Backward) {
2154 KokkosSparse::Experimental::backward_sweep_gauss_seidel_apply
2155 (mtKernelHandle_.getRawPtr(), A_->getLocalNumRows(), A_->getLocalNumCols(),
2156 kcsr.graph.row_map,kcsr.graph.entries, kcsr.values,
2157 X_colMap->getLocalViewDevice(Tpetra::Access::ReadWrite),
2158 B_in->getLocalViewDevice(Tpetra::Access::ReadOnly),
2159 zero_x_vector, update_y_vector, DampingFactor_, 1);
2162 TEUCHOS_TEST_FOR_EXCEPTION(
2163 true, std::invalid_argument,
2164 prefix <<
"The 'direction' enum does not have any of its valid "
2165 "values: Forward, Backward, or Symmetric.");
2167 update_y_vector =
false;
2170 if (copyBackOutput) {
2172 deep_copy (X , *X_domainMap);
2173 }
catch (std::exception& e) {
2174 TEUCHOS_TEST_FOR_EXCEPTION(
2175 true, std::runtime_error, prefix <<
"deep_copy(X, *X_domainMap) "
2176 "threw an exception: " << e.what ());
2180 const double dampingFlops = (DampingFactor_ == STS::one ()) ? 0.0 : 1.0;
2181 const double numVectors = as<double> (X.getNumVectors ());
2182 const double numGlobalRows = as<double> (A_->getGlobalNumRows ());
2183 const double numGlobalNonzeros = as<double> (A_->getGlobalNumEntries ());
2184 double ApplyFlops = NumSweeps_ * numVectors *
2185 (2.0 * numGlobalRows + 2.0 * numGlobalNonzeros + dampingFlops);
2186 if (direction == Tpetra::Symmetric)
2188 ApplyFlops_ += ApplyFlops;
2191template<
class MatrixType>
2194 std::ostringstream os;
2199 os <<
"\"Ifpack2::Relaxation\": {";
2201 os <<
"Initialized: " << (
isInitialized () ?
"true" :
"false") <<
", "
2202 <<
"Computed: " << (
isComputed () ?
"true" :
"false") <<
", ";
2208 if (PrecType_ == Ifpack2::Details::JACOBI) {
2210 }
else if (PrecType_ == Ifpack2::Details::GS) {
2211 os <<
"Gauss-Seidel";
2212 }
else if (PrecType_ == Ifpack2::Details::SGS) {
2213 os <<
"Symmetric Gauss-Seidel";
2214 }
else if (PrecType_ == Ifpack2::Details::MTGS) {
2215 os <<
"MT Gauss-Seidel";
2216 }
else if (PrecType_ == Ifpack2::Details::MTSGS) {
2217 os <<
"MT Symmetric Gauss-Seidel";
2218 }
else if (PrecType_ == Ifpack2::Details::GS2) {
2219 os <<
"Two-stage Gauss-Seidel";
2220 }
else if (PrecType_ == Ifpack2::Details::SGS2) {
2221 os <<
"Two-stage Symmetric Gauss-Seidel";
2226 if(hasBlockCrsMatrix_)
2229 os <<
", " <<
"sweeps: " << NumSweeps_ <<
", "
2230 <<
"damping factor: " << DampingFactor_ <<
", ";
2232 if (PrecType_ == Ifpack2::Details::MTGS || PrecType_ == Ifpack2::Details::MTSGS) {
2233 os <<
"\"relaxation: mtgs cluster size\": " << clusterSize_ <<
", ";
2234 os <<
"\"relaxation: long row threshold\": " << longRowThreshold_ <<
", ";
2235 os <<
"\"relaxation: symmetric matrix structure\": " << (is_matrix_structurally_symmetric_ ?
"true" :
"false") <<
", ";
2236 os <<
"\"relaxation: relaxation: mtgs coloring algorithm\": ";
2237 switch(mtColoringAlgorithm_)
2239 case KokkosGraph::COLORING_DEFAULT:
2240 os <<
"DEFAULT";
break;
2241 case KokkosGraph::COLORING_SERIAL:
2242 os <<
"SERIAL";
break;
2243 case KokkosGraph::COLORING_VB:
2245 case KokkosGraph::COLORING_VBBIT:
2246 os <<
"VBBIT";
break;
2247 case KokkosGraph::COLORING_VBCS:
2248 os <<
"VBCS";
break;
2249 case KokkosGraph::COLORING_VBD:
2251 case KokkosGraph::COLORING_VBDBIT:
2252 os <<
"VBDBIT";
break;
2253 case KokkosGraph::COLORING_EB:
2261 if (PrecType_ == Ifpack2::Details::GS2 ||
2262 PrecType_ == Ifpack2::Details::SGS2) {
2263 os <<
"outer sweeps: " << NumOuterSweeps_ <<
", "
2264 <<
"inner sweeps: " << NumInnerSweeps_ <<
", "
2265 <<
"inner damping factor: " << InnerDampingFactor_ <<
", ";
2269 os <<
"use l1: " << DoL1Method_ <<
", "
2270 <<
"l1 eta: " << L1Eta_ <<
", ";
2273 if (A_.is_null ()) {
2274 os <<
"Matrix: null";
2277 os <<
"Global matrix dimensions: ["
2278 << A_->getGlobalNumRows () <<
", " << A_->getGlobalNumCols () <<
"]"
2279 <<
", Global nnz: " << A_->getGlobalNumEntries();
2287template<
class MatrixType>
2290describe (Teuchos::FancyOStream &out,
2291 const Teuchos::EVerbosityLevel verbLevel)
const
2293 using Teuchos::OSTab;
2294 using Teuchos::TypeNameTraits;
2295 using Teuchos::VERB_DEFAULT;
2296 using Teuchos::VERB_NONE;
2297 using Teuchos::VERB_LOW;
2298 using Teuchos::VERB_MEDIUM;
2299 using Teuchos::VERB_HIGH;
2300 using Teuchos::VERB_EXTREME;
2303 const Teuchos::EVerbosityLevel vl =
2304 (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
2306 const int myRank = this->
getComm ()->getRank ();
2314 if (vl != VERB_NONE && myRank == 0) {
2318 out <<
"\"Ifpack2::Relaxation\":" << endl;
2320 out <<
"MatrixType: \"" << TypeNameTraits<MatrixType>::name () <<
"\""
2322 if (this->getObjectLabel () !=
"") {
2323 out <<
"Label: " << this->getObjectLabel () << endl;
2325 out <<
"Initialized: " << (
isInitialized () ?
"true" :
"false") << endl
2326 <<
"Computed: " << (
isComputed () ?
"true" :
"false") << endl
2327 <<
"Parameters: " << endl;
2330 out <<
"\"relaxation: type\": ";
2331 if (PrecType_ == Ifpack2::Details::JACOBI) {
2333 }
else if (PrecType_ == Ifpack2::Details::GS) {
2334 out <<
"Gauss-Seidel";
2335 }
else if (PrecType_ == Ifpack2::Details::SGS) {
2336 out <<
"Symmetric Gauss-Seidel";
2337 }
else if (PrecType_ == Ifpack2::Details::MTGS) {
2338 out <<
"MT Gauss-Seidel";
2339 }
else if (PrecType_ == Ifpack2::Details::MTSGS) {
2340 out <<
"MT Symmetric Gauss-Seidel";
2341 }
else if (PrecType_ == Ifpack2::Details::GS2) {
2342 out <<
"Two-stage Gauss-Seidel";
2343 }
else if (PrecType_ == Ifpack2::Details::SGS2) {
2344 out <<
"Two-stage Symmetric Gauss-Seidel";
2351 <<
"\"relaxation: sweeps\": " << NumSweeps_ << endl
2352 <<
"\"relaxation: damping factor\": " << DampingFactor_ << endl
2353 <<
"\"relaxation: min diagonal value\": " << MinDiagonalValue_ << endl
2354 <<
"\"relaxation: zero starting solution\": " << ZeroStartingSolution_ << endl
2355 <<
"\"relaxation: backward mode\": " << DoBackwardGS_ << endl
2356 <<
"\"relaxation: use l1\": " << DoL1Method_ << endl
2357 <<
"\"relaxation: l1 eta\": " << L1Eta_ << endl;
2358 if (PrecType_ == Ifpack2::Details::MTGS || PrecType_ == Ifpack2::Details::MTSGS) {
2359 out <<
"\"relaxation: mtgs cluster size\": " << clusterSize_ << endl;
2360 out <<
"\"relaxation: long row threshold\": " << longRowThreshold_ << endl;
2361 out <<
"\"relaxation: symmetric matrix structure\": " << (is_matrix_structurally_symmetric_ ?
"true" :
"false") << endl;
2362 out <<
"\"relaxation: relaxation: mtgs coloring algorithm\": ";
2363 switch(mtColoringAlgorithm_)
2365 case KokkosGraph::COLORING_DEFAULT:
2366 out <<
"DEFAULT";
break;
2367 case KokkosGraph::COLORING_SERIAL:
2368 out <<
"SERIAL";
break;
2369 case KokkosGraph::COLORING_VB:
2371 case KokkosGraph::COLORING_VBBIT:
2372 out <<
"VBBIT";
break;
2373 case KokkosGraph::COLORING_VBCS:
2374 out <<
"VBCS";
break;
2375 case KokkosGraph::COLORING_VBD:
2376 out <<
"VBD";
break;
2377 case KokkosGraph::COLORING_VBDBIT:
2378 out <<
"VBDBIT";
break;
2379 case KokkosGraph::COLORING_EB:
2386 if (PrecType_ == Ifpack2::Details::GS2 || PrecType_ == Ifpack2::Details::SGS2) {
2387 out <<
"\"relaxation: inner damping factor\": " << InnerDampingFactor_ << endl;
2388 out <<
"\"relaxation: outer sweeps\" : " << NumOuterSweeps_ << endl;
2389 out <<
"\"relaxation: inner sweeps\" : " << NumInnerSweeps_ << endl;
2392 out <<
"Computed quantities:" << endl;
2395 out <<
"Global number of rows: " << A_->getGlobalNumRows () << endl
2396 <<
"Global number of columns: " << A_->getGlobalNumCols () << endl;
2399 out <<
"Properties of input diagonal entries:" << endl;
2402 out <<
"Magnitude of minimum-magnitude diagonal entry: "
2403 << globalMinMagDiagEntryMag_ << endl
2404 <<
"Magnitude of maximum-magnitude diagonal entry: "
2405 << globalMaxMagDiagEntryMag_ << endl
2406 <<
"Number of diagonal entries with small magnitude: "
2407 << globalNumSmallDiagEntries_ << endl
2408 <<
"Number of zero diagonal entries: "
2409 << globalNumZeroDiagEntries_ << endl
2410 <<
"Number of diagonal entries with negative real part: "
2411 << globalNumNegDiagEntries_ << endl
2412 <<
"Abs 2-norm diff between computed and actual inverse "
2413 <<
"diagonal: " << globalDiagNormDiff_ << endl;
2417 out <<
"Saved diagonal offsets: "
2418 << (savedDiagOffsets_ ?
"true" :
"false") << endl;
2420 out <<
"Call counts and total times (in seconds): " << endl;
2423 out <<
"initialize: " << endl;
2426 out <<
"Call count: " << NumInitialize_ << endl;
2427 out <<
"Total time: " << InitializeTime_ << endl;
2429 out <<
"compute: " << endl;
2432 out <<
"Call count: " << NumCompute_ << endl;
2433 out <<
"Total time: " << ComputeTime_ << endl;
2435 out <<
"apply: " << endl;
2438 out <<
"Call count: " << NumApply_ << endl;
2439 out <<
"Total time: " << ApplyTime_ << endl;
2448#define IFPACK2_RELAXATION_INSTANT(S,LO,GO,N) \
2449 template class Ifpack2::Relaxation< Tpetra::RowMatrix<S, LO, GO, N> >;
File for utility functions.
Compute scaled damped residual for Chebyshev.
Definition Ifpack2_Details_InverseDiagonalKernel_decl.hpp:77
Diagonal preconditioner.
Definition Ifpack2_Diagonal_decl.hpp:80
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition Ifpack2_Relaxation_decl.hpp:248
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition Ifpack2_Relaxation_def.hpp:214
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_Relaxation_decl.hpp:260
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition Ifpack2_Relaxation_def.hpp:467
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_Relaxation_decl.hpp:263
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition Ifpack2_Relaxation_def.hpp:477
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_Relaxation_def.hpp:497
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition Ifpack2_Relaxation_def.hpp:556
int getNumCompute() const
Total number of calls to compute().
Definition Ifpack2_Relaxation_def.hpp:532
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition Ifpack2_Relaxation_def.hpp:488
void compute()
Compute the preconditioner ("numeric setup");.
Definition Ifpack2_Relaxation_def.hpp:1020
void initialize()
Initialize the preconditioner ("symbolic setup").
Definition Ifpack2_Relaxation_def.hpp:722
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition Ifpack2_Relaxation_def.hpp:550
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object's attributes to the given output stream.
Definition Ifpack2_Relaxation_def.hpp:2290
std::string description() const
A simple one-line description of this object.
Definition Ifpack2_Relaxation_def.hpp:2192
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition Ifpack2_Relaxation_def.hpp:226
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition Ifpack2_Relaxation_def.hpp:544
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_Relaxation_def.hpp:193
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_Relaxation_decl.hpp:254
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 preconditioner to X, returning the result in Y.
Definition Ifpack2_Relaxation_def.hpp:588
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_Relaxation_decl.hpp:257
bool isComputed() const
Return true if compute() has been called.
Definition Ifpack2_Relaxation_decl.hpp:426
void applyMat(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) const
Apply the input matrix to X, returning the result in Y.
Definition Ifpack2_Relaxation_def.hpp:701
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition Ifpack2_Relaxation_decl.hpp:273
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_Relaxation_decl.hpp:269
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_Relaxation_def.hpp:575
int getNumApply() const
Total number of calls to apply().
Definition Ifpack2_Relaxation_def.hpp:538
int getNumInitialize() const
Total number of calls to initialize().
Definition Ifpack2_Relaxation_def.hpp:526
double getComputeFlops() const
Total number of floating-point operations over all calls to compute().
Definition Ifpack2_Relaxation_def.hpp:562
bool hasTransposeApply() const
Whether apply() and applyMat() let you apply the transpose or conjugate transpose.
Definition Ifpack2_Relaxation_def.hpp:520
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition Ifpack2_Relaxation_decl.hpp:418
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_Relaxation_def.hpp:510
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition Ifpack2_Relaxation_def.hpp:568
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74
void getValidParameters(Teuchos::ParameterList ¶ms)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition Ifpack2_Parameters.cpp:51