43 #ifndef IFPACK2_BLOCKRELAXATION_DEF_HPP
44 #define IFPACK2_BLOCKRELAXATION_DEF_HPP
47 #include "Ifpack2_LinearPartitioner.hpp"
48 #include "Ifpack2_LinePartitioner.hpp"
50 #include "Ifpack2_Details_UserPartitioner_def.hpp"
51 #include <Ifpack2_Parameters.hpp>
55 template<
class MatrixType,
class ContainerType>
60 IsInitialized_ =
false;
62 Partitioner_ = Teuchos::null;
66 IsParallel_ = (A->getRowMap ()->getComm ()->getSize () != 1);
68 Teuchos::rcp_dynamic_cast<
const block_crs_matrix_type> (A);
69 hasBlockCrsMatrix_ = !A_bcrs.
is_null();
71 if (! Container_.is_null ()) {
72 Container_->clearBlocks();
81 template<
class MatrixType,
class ContainerType>
85 Time_ (Teuchos::
rcp (new Teuchos::Time (
"Ifpack2::BlockRelaxation"))),
86 Container_ (Teuchos::null),
87 Partitioner_ (Teuchos::null),
88 PartitionerType_ (
"linear"),
91 containerType_ (
"TriDi"),
92 PrecType_ (Ifpack2::Details::JACOBI),
93 ZeroStartingSolution_ (true),
94 hasBlockCrsMatrix_ (false),
95 DoBackwardGS_ (false),
97 DampingFactor_ (
STS::one ()),
98 IsInitialized_ (false),
103 InitializeTime_ (0.0),
107 NumGlobalNonzeros_ (0),
109 Importer_ (Teuchos::null)
114 template<
class MatrixType,
class ContainerType>
119 template<
class MatrixType,
class ContainerType>
126 validParams->
set(
"relaxation: container",
"TriDi");
127 validParams->
set(
"relaxation: backward mode",
false);
128 validParams->
set(
"relaxation: type",
"Jacobi");
129 validParams->
set(
"relaxation: sweeps", (
int)1);
130 validParams->
set(
"relaxation: damping factor", STS::one());
131 validParams->
set(
"relaxation: zero starting solution",
true);
132 validParams->
set(
"schwarz: compute condest",
false);
133 validParams->
set(
"schwarz: combine mode",
"ZERO");
134 validParams->
set(
"schwarz: use reordering",
true);
135 validParams->
set(
"schwarz: filter singletons",
false);
136 validParams->
set(
"schwarz: overlap level", (
int)0);
137 validParams->
set(
"partitioner: type",
"greedy");
138 validParams->
set(
"partitioner: local parts", (
int)1);
139 validParams->
set(
"partitioner: overlap", (
int)0);
141 validParams->
set(
"partitioner: parts", tmp0);
142 validParams->
set(
"partitioner: maintain sparsity",
false);
145 validParams->
set(
"Amesos2",dummyList);
147 validParams->
set(
"Amesos2 solver name",
"KLU2");
150 validParams->
set(
"partitioner: map", tmp);
152 validParams->
set(
"partitioner: line detection threshold",(
double)0.0);
153 validParams->
set(
"partitioner: PDE equations",(
int)1);
155 typename MatrixType::local_ordinal_type,
156 typename MatrixType::global_ordinal_type,
157 typename MatrixType::node_type> > dummy;
158 validParams->
set(
"partitioner: coordinates",dummy);
163 template<
class MatrixType,
class ContainerType>
179 containerType_ = List.
get<std::string> (
"relaxation: container");
181 if(containerType_ ==
"Tridiagonal") {
182 containerType_ =
"TriDi";
184 if(containerType_ ==
"Block Tridiagonal") {
185 containerType_ =
"BlockTriDi";
190 const std::string relaxType = List.
get<std::string> (
"relaxation: type");
192 if (relaxType ==
"Jacobi") {
193 PrecType_ = Ifpack2::Details::JACOBI;
195 else if (relaxType ==
"MT Split Jacobi") {
196 PrecType_ = Ifpack2::Details::MTSPLITJACOBI;
198 else if (relaxType ==
"Gauss-Seidel") {
199 PrecType_ = Ifpack2::Details::GS;
201 else if (relaxType ==
"Symmetric Gauss-Seidel") {
202 PrecType_ = Ifpack2::Details::SGS;
206 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
207 "setParameters: Invalid parameter value \"" << relaxType
208 <<
"\" for parameter \"relaxation: type\".");
213 NumSweeps_ = List.
get<
int> (
"relaxation: sweeps");
219 if (List.
isParameter (
"relaxation: damping factor")) {
220 if (List.
isType<
double> (
"relaxation: damping factor")) {
221 const double dampingFactor =
222 List.
get<
double> (
"relaxation: damping factor");
223 DampingFactor_ =
static_cast<scalar_type> (dampingFactor);
226 DampingFactor_ = List.
get<
scalar_type> (
"relaxation: damping factor");
231 DampingFactor_ =
static_cast<scalar_type> (dampingFactor);
235 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
236 "setParameters: Parameter \"relaxation: damping factor\" "
237 "has the wrong type.");
241 if (List.
isParameter (
"relaxation: zero starting solution")) {
242 ZeroStartingSolution_ = List.
get<
bool> (
"relaxation: zero starting solution");
245 if (List.
isParameter (
"relaxation: backward mode")) {
246 DoBackwardGS_ = List.
get<
bool> (
"relaxation: backward mode");
250 PartitionerType_ = List.
get<std::string> (
"partitioner: type");
255 if (List.
isParameter (
"partitioner: local parts")) {
259 else if (! std::is_same<int, local_ordinal_type>::value &&
260 List.
isType<
int> (
"partitioner: local parts")) {
261 NumLocalBlocks_ = List.
get<
int> (
"partitioner: local parts");
265 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
266 "setParameters: Parameter \"partitioner: local parts\" "
267 "has the wrong type.");
271 if (List.
isParameter (
"partitioner: overlap level")) {
272 if (List.
isType<
int> (
"partitioner: overlap level")) {
273 OverlapLevel_ = List.
get<
int> (
"partitioner: overlap level");
275 else if (! std::is_same<int, local_ordinal_type>::value &&
281 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
282 "setParameters: Parameter \"partitioner: overlap level\" "
283 "has the wrong type.");
287 std::string defaultContainer =
"TriDi";
294 if (PrecType_ != Ifpack2::Details::JACOBI) {
297 if (NumLocalBlocks_ < static_cast<local_ordinal_type> (0)) {
298 NumLocalBlocks_ = A_->getNodeNumRows() / (-NumLocalBlocks_);
304 DoBackwardGS_, std::runtime_error,
305 "Ifpack2::BlockRelaxation:setParameters: Setting the \"relaxation: "
306 "backward mode\" parameter to true is not yet supported.");
313 template<
class MatrixType,
class ContainerType>
318 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::getComm: "
319 "The matrix is null. You must call setMatrix() with a nonnull matrix "
320 "before you may call this method.");
321 return A_->getComm ();
324 template<
class MatrixType,
class ContainerType>
325 Teuchos::RCP<
const Tpetra::RowMatrix<
typename MatrixType::scalar_type,
326 typename MatrixType::local_ordinal_type,
327 typename MatrixType::global_ordinal_type,
328 typename MatrixType::node_type> >
333 template<
class MatrixType,
class ContainerType>
334 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
335 typename MatrixType::global_ordinal_type,
336 typename MatrixType::node_type> >
341 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
342 "getDomainMap: The matrix is null. You must call setMatrix() with a "
343 "nonnull matrix before you may call this method.");
344 return A_->getDomainMap ();
347 template<
class MatrixType,
class ContainerType>
348 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
349 typename MatrixType::global_ordinal_type,
350 typename MatrixType::node_type> >
355 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
356 "getRangeMap: The matrix is null. You must call setMatrix() with a "
357 "nonnull matrix before you may call this method.");
358 return A_->getRangeMap ();
361 template<
class MatrixType,
class ContainerType>
368 template<
class MatrixType,
class ContainerType>
372 return NumInitialize_;
375 template<
class MatrixType,
class ContainerType>
383 template<
class MatrixType,
class ContainerType>
391 template<
class MatrixType,
class ContainerType>
396 return InitializeTime_;
399 template<
class MatrixType,
class ContainerType>
407 template<
class MatrixType,
class ContainerType>
416 template<
class MatrixType,
class ContainerType>
419 A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::getNodeSmootherComplexity: "
420 "The input matrix A is null. Please call setMatrix() with a nonnull "
421 "input matrix, then call compute(), before calling this method.");
424 size_t block_nnz = 0;
426 block_nnz += Partitioner_->numRowsInPart(i) *Partitioner_->numRowsInPart(i);
428 return block_nnz + A_->getNodeNumEntries();
431 template<
class MatrixType,
class ContainerType>
434 apply (
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
435 typename MatrixType::local_ordinal_type,
436 typename MatrixType::global_ordinal_type,
437 typename MatrixType::node_type>& X,
438 Tpetra::MultiVector<
typename MatrixType::scalar_type,
439 typename MatrixType::local_ordinal_type,
440 typename MatrixType::global_ordinal_type,
441 typename MatrixType::node_type>& Y,
447 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::apply: "
448 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
449 "then call initialize() and compute() (in that order), before you may "
450 "call this method.");
452 ! isComputed (), std::runtime_error,
"Ifpack2::BlockRelaxation::apply: "
453 "isComputed() must be true prior to calling apply.");
454 TEUCHOS_TEST_FOR_EXCEPTION(
455 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
456 "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = "
457 << X.getNumVectors () <<
" != Y.getNumVectors() = "
458 << Y.getNumVectors () <<
".");
459 TEUCHOS_TEST_FOR_EXCEPTION(
461 "apply: This method currently only implements the case mode == Teuchos::"
462 "NO_TRANS. Transposed modes are not currently supported.");
463 TEUCHOS_TEST_FOR_EXCEPTION(
465 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
466 "the case alpha == 1. You specified alpha = " << alpha <<
".");
467 TEUCHOS_TEST_FOR_EXCEPTION(
469 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
470 "the case beta == 0. You specified beta = " << beta <<
".");
478 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
479 auto X_lcl_host = X.template getLocalView<Kokkos::HostSpace> ();
480 auto Y_lcl_host = Y.template getLocalView<Kokkos::HostSpace> ();
482 auto X_lcl_host = X.getLocalViewHost ();
483 auto Y_lcl_host = Y.getLocalViewHost ();
485 if (X_lcl_host.data () == Y_lcl_host.data ()) {
488 X_copy = rcpFromRef (X);
492 if (ZeroStartingSolution_) {
493 Y.putScalar (STS::zero ());
498 case Ifpack2::Details::JACOBI:
499 ApplyInverseJacobi(*X_copy,Y);
501 case Ifpack2::Details::GS:
502 ApplyInverseGS(*X_copy,Y);
504 case Ifpack2::Details::SGS:
505 ApplyInverseSGS(*X_copy,Y);
507 case Ifpack2::Details::MTSPLITJACOBI:
508 Container_->applyInverseJacobi(*X_copy, Y, ZeroStartingSolution_, NumSweeps_);
511 TEUCHOS_TEST_FOR_EXCEPTION
512 (
true, std::logic_error,
"Ifpack2::BlockRelaxation::apply: Invalid "
513 "PrecType_ enum value " << PrecType_ <<
". Valid values are Ifpack2::"
514 "Details::JACOBI = " << Ifpack2::Details::JACOBI <<
", Ifpack2::Details"
515 "::GS = " << Ifpack2::Details::GS <<
", and Ifpack2::Details::SGS = "
516 << Ifpack2::Details::SGS <<
". Please report this bug to the Ifpack2 "
522 ApplyTime_ += Time_->totalElapsedTime();
525 template<
class MatrixType,
class ContainerType>
528 applyMat (
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
529 typename MatrixType::local_ordinal_type,
530 typename MatrixType::global_ordinal_type,
531 typename MatrixType::node_type>& X,
532 Tpetra::MultiVector<
typename MatrixType::scalar_type,
533 typename MatrixType::local_ordinal_type,
534 typename MatrixType::global_ordinal_type,
535 typename MatrixType::node_type>& Y,
538 A_->apply (X, Y, mode);
541 template<
class MatrixType,
class ContainerType>
547 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
551 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::initialize: "
552 "The matrix is null. You must call setMatrix() with a nonnull matrix "
553 "before you may call this method.");
557 Teuchos::rcp_dynamic_cast<
const block_crs_matrix_type> (A_);
558 hasBlockCrsMatrix_ = !A_bcrs.
is_null();
559 if (A_bcrs.is_null ()) {
560 hasBlockCrsMatrix_ =
false;
563 hasBlockCrsMatrix_ =
true;
566 IsInitialized_ =
false;
569 NumLocalRows_ = A_->getNodeNumRows ();
570 NumGlobalRows_ = A_->getGlobalNumRows ();
571 NumGlobalNonzeros_ = A_->getGlobalNumEntries ();
576 Partitioner_ = Teuchos::null;
578 if (PartitionerType_ ==
"linear") {
581 }
else if (PartitionerType_ ==
"line") {
584 }
else if (PartitionerType_ ==
"user") {
591 (
true, std::logic_error,
"Ifpack2::BlockRelaxation::initialize: Unknown "
592 "partitioner type " << PartitionerType_ <<
". Valid values are "
593 "\"linear\", \"line\", and \"user\".");
597 Partitioner_->setParameters (List_);
598 Partitioner_->compute ();
601 NumLocalBlocks_ = Partitioner_->numLocalParts ();
606 if (A_->getComm()->getSize() != 1) {
615 (NumSweeps_ < 0, std::logic_error,
"Ifpack2::BlockRelaxation::initialize: "
616 "NumSweeps_ = " << NumSweeps_ <<
" < 0.");
619 ExtractSubmatricesStructure ();
623 if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) {
625 (hasBlockCrsMatrix_, std::runtime_error,
626 "Ifpack2::BlockRelaxation::initialize: "
627 "We do not support overlapped Jacobi yet for Tpetra::BlockCrsMatrix. Sorry!");
630 W_ =
rcp (
new vector_type (A_->getRowMap ()));
631 W_->putScalar (STS::zero ());
635 for (
size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) {
638 int LID = (*Partitioner_)(i,j);
639 w_ptr[LID]+= STS::one();
642 W_->reciprocal (*W_);
647 InitializeTime_ += Time_->totalElapsedTime ();
648 IsInitialized_ =
true;
652 template<
class MatrixType,
class ContainerType>
660 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::compute: "
661 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
662 "then call initialize(), before you may call this method.");
664 if (! isInitialized ()) {
673 Container_->compute();
677 ComputeTime_ += Time_->totalElapsedTime();
681 template<
class MatrixType,
class ContainerType>
687 (Partitioner_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
688 "ExtractSubmatricesStructure: Partitioner object is null.");
690 NumLocalBlocks_ = Partitioner_->numLocalParts ();
691 std::string containerType = ContainerType::getName ();
692 if (containerType ==
"Generic") {
696 containerType = containerType_;
700 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) {
701 const size_t numRows = Partitioner_->numRowsInPart (i);
703 localRows[i].resize(numRows);
705 for (
size_t j = 0; j < numRows; ++j) {
706 localRows[i][j] = (*Partitioner_) (i,j);
710 Container_->setParameters(List_);
711 Container_->initialize();
714 template<
class MatrixType,
class ContainerType>
716 BlockRelaxation<MatrixType,ContainerType>::
717 ApplyInverseJacobi (
const MV& X, MV& Y)
const
719 const size_t NumVectors = X.getNumVectors ();
720 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
721 typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace> ();
722 typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace> ();
724 auto XView = X.getLocalViewHost ();
725 auto YView = Y.getLocalViewHost ();
727 MV AY (Y.getMap (), NumVectors);
729 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
730 typename ContainerType::HostView AYView = AY.template getLocalView<Kokkos::HostSpace> ();
732 auto AYView = AY.getLocalViewHost ();
735 int starting_iteration = 0;
736 if (OverlapLevel_ > 0)
739 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
740 typename ContainerType::HostView WView = W_->template getLocalView<Kokkos::HostSpace> ();
742 auto WView = W_->getLocalViewHost ();
744 if(ZeroStartingSolution_) {
745 Container_->DoOverlappingJacobi(XView, YView, WView, X.getStride());
746 starting_iteration = 1;
748 const scalar_type ONE = STS::one();
749 for(
int j = starting_iteration; j < NumSweeps_; j++)
752 AY.update (ONE, X, -ONE);
753 Container_->DoOverlappingJacobi (AYView, YView, WView, X.getStride());
759 if(ZeroStartingSolution_)
761 Container_->DoJacobi (XView, YView, X.getStride());
762 starting_iteration = 1;
764 const scalar_type ONE = STS::one();
765 for(
int j = starting_iteration; j < NumSweeps_; j++)
768 AY.update (ONE, X, -ONE);
769 Container_->DoJacobi (AYView, YView, X.getStride());
774 template<
class MatrixType,
class ContainerType>
776 BlockRelaxation<MatrixType,ContainerType>::
777 ApplyInverseGS (
const MV& X, MV& Y)
const
781 size_t numVecs = X.getNumVectors();
783 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
784 typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace> ();
785 typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace> ();
787 auto XView = X.getLocalViewHost ();
788 auto YView = Y.getLocalViewHost ();
792 bool deleteY2 =
false;
795 Y2 = ptr(
new MV(Importer_->getTargetMap(), numVecs));
802 for(
int j = 0; j < NumSweeps_; ++j)
805 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
806 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
807 typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace> ();
809 auto Y2View = Y2->getLocalViewHost ();
811 Container_->DoGaussSeidel(XView, YView, Y2View, X.getStride());
816 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
817 typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace> ();
819 auto Y2View = Y2->getLocalViewHost ();
821 for(
int j = 0; j < NumSweeps_; ++j)
823 Container_->DoGaussSeidel(XView, YView, Y2View, X.getStride());
830 template<
class MatrixType,
class ContainerType>
832 BlockRelaxation<MatrixType,ContainerType>::
833 ApplyInverseSGS (
const MV& X, MV& Y)
const
838 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
839 typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace> ();
840 typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace> ();
842 auto XView = X.getLocalViewHost ();
843 auto YView = Y.getLocalViewHost ();
847 bool deleteY2 =
false;
850 Y2 = ptr(
new MV(Importer_->getTargetMap(), X.getNumVectors()));
857 for(
int j = 0; j < NumSweeps_; ++j)
860 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
861 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
862 typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace> ();
864 auto Y2View = Y2->getLocalViewHost ();
866 Container_->DoSGS(XView, YView, Y2View, X.getStride());
871 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
872 typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace> ();
874 auto Y2View = Y2->getLocalViewHost ();
876 for(
int j = 0; j < NumSweeps_; ++j)
878 Container_->DoSGS(XView, YView, Y2View, X.getStride());
885 template<
class MatrixType,
class ContainerType>
886 void BlockRelaxation<MatrixType,ContainerType>::computeImporter ()
const
894 using Teuchos::rcp_dynamic_cast;
897 if(hasBlockCrsMatrix_)
899 const RCP<const block_crs_matrix_type> bcrs =
900 rcp_dynamic_cast<
const block_crs_matrix_type>(A_);
901 int bs_ = bcrs->getBlockSize();
902 RCP<const map_type> oldDomainMap = A_->getDomainMap();
903 RCP<const map_type> oldColMap = A_->getColMap();
906 global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
907 global_ordinal_type indexBase = oldColMap->getIndexBase();
908 RCP<const Comm<int>> comm = oldColMap->getComm();
909 ArrayView<const global_ordinal_type> oldColElements = oldColMap->getNodeElementList();
910 Array<global_ordinal_type> newColElements(bs_ * oldColElements.size());
911 for(
int i = 0; i < oldColElements.size(); i++)
913 for(
int j = 0; j < bs_; j++)
914 newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
916 RCP<map_type> colMap(
new map_type(numGlobalElements, newColElements, indexBase, comm));
918 Importer_ =
rcp(
new import_type(oldDomainMap, colMap));
920 else if(!A_.is_null())
922 Importer_ = A_->getGraph()->getImporter();
923 if(Importer_.is_null())
924 Importer_ =
rcp(
new import_type(A_->getDomainMap(), A_->getColMap()));
930 template<
class MatrixType,
class ContainerType>
935 std::ostringstream out;
940 out <<
"\"Ifpack2::BlockRelaxation\": {";
941 if (this->getObjectLabel () !=
"") {
942 out <<
"Label: \"" << this->getObjectLabel () <<
"\", ";
947 out <<
"Matrix: null, ";
952 out <<
"Global matrix dimensions: ["
953 << A_->getGlobalNumRows () <<
", " << A_->getGlobalNumCols () <<
"], ";
958 out <<
"\"relaxation: type\": ";
959 if (PrecType_ == Ifpack2::Details::JACOBI) {
960 out <<
"Block Jacobi";
961 }
else if (PrecType_ == Ifpack2::Details::GS) {
962 out <<
"Block Gauss-Seidel";
963 }
else if (PrecType_ == Ifpack2::Details::SGS) {
964 out <<
"Block Symmetric Gauss-Seidel";
969 out <<
", overlap: " << OverlapLevel_;
971 out <<
", " <<
"sweeps: " << NumSweeps_ <<
", "
972 <<
"damping factor: " << DampingFactor_ <<
", ";
974 std::string containerType = ContainerType::getName();
975 out <<
"block container: " << ((containerType ==
"Generic") ? containerType_ : containerType);
976 if(List_.isParameter(
"partitioner: PDE equations"))
977 out <<
", dofs/node: "<<List_.get<
int>(
"partitioner: PDE equations");
984 template<
class MatrixType,
class ContainerType>
1002 const int myImageID = A_->getComm()->getRank();
1012 if (vl !=
VERB_NONE && myImageID == 0) {
1013 out <<
"Ifpack2::BlockRelaxation<"
1014 << TypeNameTraits<MatrixType>::name () <<
", "
1015 << TypeNameTraits<ContainerType>::name () <<
" >:";
1018 if (this->getObjectLabel () !=
"") {
1019 out <<
"label: \"" << this->getObjectLabel () <<
"\"" << endl;
1021 out <<
"initialized: " << (isInitialized () ?
"true" :
"false") << endl
1022 <<
"computed: " << (isComputed () ?
"true" :
"false") << endl;
1024 std::string precType;
1025 if (PrecType_ == Ifpack2::Details::JACOBI) {
1026 precType =
"Block Jacobi";
1027 }
else if (PrecType_ == Ifpack2::Details::GS) {
1028 precType =
"Block Gauss-Seidel";
1029 }
else if (PrecType_ == Ifpack2::Details::SGS) {
1030 precType =
"Block Symmetric Gauss-Seidel";
1032 out <<
"type: " << precType << endl
1033 <<
"global number of rows: " << A_->getGlobalNumRows () << endl
1034 <<
"global number of columns: " << A_->getGlobalNumCols () << endl
1035 <<
"number of sweeps: " << NumSweeps_ << endl
1036 <<
"damping factor: " << DampingFactor_ << endl
1037 <<
"backwards mode: "
1038 << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ?
"true" :
"false")
1040 <<
"zero starting solution: "
1041 << (ZeroStartingSolution_ ?
"true" :
"false") << endl;
1059 #ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1061 #define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \
1063 class Ifpack2::BlockRelaxation< \
1064 Tpetra::RowMatrix<S, LO, GO, N>, \
1065 Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N> > >;
1067 #endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1069 #endif // IFPACK2_BLOCKRELAXATION_DEF_HPP
Ifpack2::BlockRelaxation class declaration.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition: Ifpack2_BlockRelaxation_def.hpp:315
ParameterList & disableRecursiveValidation()
std::string description() const
A one-line description of this object.
Definition: Ifpack2_BlockRelaxation_def.hpp:933
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:410
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static Teuchos::RCP< BaseContainer > build(std::string containerType, const Teuchos::RCP< const MatrixType > &A, const Teuchos::Array< Teuchos::Array< local_ordinal_type >> &localRows, const Teuchos::RCP< const import_type > importer, int OverlapLevel, scalar_type DampingFactor)
Build a specialization of Ifpack2::Container given a key that has been registered.
Definition: Ifpack2_ContainerFactory_def.hpp:89
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:386
void getParameter(const Teuchos::ParameterList ¶ms, const std::string &name, T &value)
Set a value from a ParameterList if a parameter with the specified name exists.
Definition: Ifpack2_Parameters.hpp:59
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:352
int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:378
void setParameters(const Teuchos::ParameterList ¶ms)
Sets all the parameters for the preconditioner.
Definition: Ifpack2_BlockRelaxation_def.hpp:166
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition: Ifpack2_BlockRelaxation_def.hpp:655
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition: Ifpack2_BlockRelaxation_decl.hpp:83
double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:394
Declaration of a user-defined partitioner in which the user can define a partition of the graph...
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
double getComputeTime() const
Returns the time spent in compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:402
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_BlockRelaxation_def.hpp:417
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_BlockRelaxation_decl.hpp:109
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_BlockRelaxation_def.hpp:57
virtual ~BlockRelaxation()
Destructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:116
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:338
Partition in which the user can define a nonoverlapping partition of the graph in any way they choose...
Definition: Ifpack2_Details_UserPartitioner_decl.hpp:69
int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:371
void initialize()
Initialize.
Definition: Ifpack2_BlockRelaxation_def.hpp:544
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
Definition: Ifpack2_BlockRelaxation_def.hpp:987
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:101
void applyMat(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS) const
Applies the matrix to a Tpetra::MultiVector.
Definition: Ifpack2_BlockRelaxation_def.hpp:528
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:98
Interface for creating and solving a local linear problem.
Definition: Ifpack2_Container.hpp:114
bool isType(const std::string &name) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner's constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:329
void getValidParameters(Teuchos::ParameterList ¶ms)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition: Ifpack2_Parameters.cpp:50
Ifpack2::LinePartitioner: A class to define partitions into a set of lines.
Definition: Ifpack2_LinePartitioner_decl.hpp:77
A class to define linear partitions.
Definition: Ifpack2_LinearPartitioner_decl.hpp:60
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_BlockRelaxation_def.hpp:122
void apply(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Applies the preconditioner to X, returns the result in Y.
Definition: Ifpack2_BlockRelaxation_def.hpp:434
BlockRelaxation(const Teuchos::RCP< const row_matrix_type > &Matrix)
Constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:83