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>
52#include "Teuchos_TimeMonitor.hpp"
56template<
class MatrixType,
class ContainerType>
58setMatrix (
const Teuchos::RCP<const row_matrix_type>& A)
60 if (A.getRawPtr () != A_.getRawPtr ()) {
61 IsInitialized_ =
false;
63 Partitioner_ = Teuchos::null;
67 IsParallel_ = (A->getRowMap ()->
getComm ()->getSize () != 1);
68 Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
69 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A);
70 hasBlockCrsMatrix_ = !A_bcrs.is_null();
72 if (! Container_.is_null ()) {
75 Container_->clearBlocks();
84template<
class MatrixType,
class ContainerType>
88 Container_ (Teuchos::null),
89 Partitioner_ (Teuchos::null),
90 PartitionerType_ (
"linear"),
93 containerType_ (
"TriDi"),
95 ZeroStartingSolution_ (true),
96 hasBlockCrsMatrix_ (false),
97 DoBackwardGS_ (false),
100 schwarzCombineMode_(
"ZERO"),
101 DampingFactor_ (STS::one ()),
102 IsInitialized_ (false),
106 TimerForApply_(true),
108 InitializeTime_ (0.0),
112 NumGlobalNonzeros_ (0),
114 Importer_ (Teuchos::null)
119template<
class MatrixType,
class ContainerType>
124template<
class MatrixType,
class ContainerType>
125Teuchos::RCP<const Teuchos::ParameterList>
129 Teuchos::RCP<Teuchos::ParameterList> validParams = Teuchos::parameterList (
"Ifpack2::BlockRelaxation");
131 validParams->set(
"relaxation: container",
"TriDi");
132 validParams->set(
"relaxation: backward mode",
false);
133 validParams->set(
"relaxation: type",
"Jacobi");
134 validParams->set(
"relaxation: sweeps", 1);
135 validParams->set(
"relaxation: damping factor", STS::one());
136 validParams->set(
"relaxation: zero starting solution",
true);
137 validParams->set(
"block relaxation: decouple dofs",
false);
138 validParams->set(
"schwarz: compute condest",
false);
139 validParams->set(
"schwarz: combine mode",
"ZERO");
140 validParams->set(
"schwarz: use reordering",
true);
141 validParams->set(
"schwarz: filter singletons",
false);
142 validParams->set(
"schwarz: overlap level", 0);
143 validParams->set(
"partitioner: type",
"greedy");
144 validParams->set(
"partitioner: local parts", 1);
145 validParams->set(
"partitioner: overlap", 0);
146 validParams->set(
"partitioner: combine mode",
"ZERO");
147 Teuchos::Array<Teuchos::ArrayRCP<int>> tmp0;
148 validParams->set(
"partitioner: parts", tmp0);
149 Teuchos::Array<Teuchos::ArrayRCP<typename MatrixType::global_ordinal_type> > tmp1;
150 validParams->set(
"partitioner: global ID parts", tmp1);
151 validParams->set(
"partitioner: nonsymmetric overlap combine",
false);
152 validParams->set(
"partitioner: maintain sparsity",
false);
153 validParams->set(
"fact: ilut level-of-fill", 1.0);
154 validParams->set(
"fact: absolute threshold", 0.0);
155 validParams->set(
"fact: relative threshold", 1.0);
156 validParams->set(
"fact: relax value", 0.0);
158 Teuchos::ParameterList dummyList;
159 validParams->set(
"Amesos2",dummyList);
160 validParams->sublist(
"Amesos2").disableRecursiveValidation();
161 validParams->set(
"Amesos2 solver name",
"KLU2");
163 Teuchos::ArrayRCP<int> tmp;
164 validParams->set(
"partitioner: map", tmp);
166 validParams->set(
"partitioner: line detection threshold", 0.0);
167 validParams->set(
"partitioner: PDE equations", 1);
168 Teuchos::RCP<Tpetra::MultiVector<typename Teuchos::ScalarTraits<typename MatrixType::scalar_type>::magnitudeType,
169 typename MatrixType::local_ordinal_type,
170 typename MatrixType::global_ordinal_type,
171 typename MatrixType::node_type> > dummy;
172 validParams->set(
"partitioner: coordinates",dummy);
173 validParams->set(
"timer for apply",
true);
178template<
class MatrixType,
class ContainerType>
186 this->setParametersImpl(
const_cast<Teuchos::ParameterList&
>(pl));
189template<
class MatrixType,
class ContainerType>
194 if (List.isType<
double>(
"relaxation: damping factor")) {
197 scalar_type df = List.get<
double>(
"relaxation: damping factor");
198 List.remove(
"relaxation: damping factor");
199 List.set(
"relaxation: damping factor",df);
203 Teuchos::RCP<const Teuchos::ParameterList> validparams;
205 List.validateParameters (*validparams);
207 if (List.isParameter (
"relaxation: container")) {
213 containerType_ = List.get<std::string> (
"relaxation: container");
215 if(containerType_ ==
"Tridiagonal") {
216 containerType_ =
"TriDi";
218 if(containerType_ ==
"Block Tridiagonal") {
219 containerType_ =
"BlockTriDi";
223 if (List.isParameter (
"relaxation: type")) {
224 const std::string relaxType = List.get<std::string> (
"relaxation: type");
226 if (relaxType ==
"Jacobi") {
227 PrecType_ = Ifpack2::Details::JACOBI;
229 else if (relaxType ==
"MT Split Jacobi") {
230 PrecType_ = Ifpack2::Details::MTSPLITJACOBI;
232 else if (relaxType ==
"Gauss-Seidel") {
233 PrecType_ = Ifpack2::Details::GS;
235 else if (relaxType ==
"Symmetric Gauss-Seidel") {
236 PrecType_ = Ifpack2::Details::SGS;
239 TEUCHOS_TEST_FOR_EXCEPTION
240 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
241 "setParameters: Invalid parameter value \"" << relaxType
242 <<
"\" for parameter \"relaxation: type\".");
246 if (List.isParameter (
"relaxation: sweeps")) {
247 NumSweeps_ = List.get<
int> (
"relaxation: sweeps");
253 if (List.isParameter (
"relaxation: damping factor")) {
254 if (List.isType<
double> (
"relaxation: damping factor")) {
255 const double dampingFactor =
256 List.get<
double> (
"relaxation: damping factor");
257 DampingFactor_ =
static_cast<scalar_type
> (dampingFactor);
259 else if (List.isType<scalar_type> (
"relaxation: damping factor")) {
260 DampingFactor_ = List.get<scalar_type> (
"relaxation: damping factor");
262 else if (List.isType<magnitude_type> (
"relaxation: damping factor")) {
263 const magnitude_type dampingFactor =
264 List.get<magnitude_type> (
"relaxation: damping factor");
265 DampingFactor_ =
static_cast<scalar_type
> (dampingFactor);
268 TEUCHOS_TEST_FOR_EXCEPTION
269 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
270 "setParameters: Parameter \"relaxation: damping factor\" "
271 "has the wrong type.");
275 if (List.isParameter (
"relaxation: zero starting solution")) {
276 ZeroStartingSolution_ = List.get<
bool> (
"relaxation: zero starting solution");
279 if (List.isParameter (
"relaxation: backward mode")) {
280 DoBackwardGS_ = List.get<
bool> (
"relaxation: backward mode");
283 if (List.isParameter (
"partitioner: type")) {
284 PartitionerType_ = List.get<std::string> (
"partitioner: type");
289 if (List.isParameter (
"partitioner: local parts")) {
293 else if (! std::is_same<int, local_ordinal_type>::value &&
294 List.isType<
int> (
"partitioner: local parts")) {
295 NumLocalBlocks_ = List.get<
int> (
"partitioner: local parts");
298 TEUCHOS_TEST_FOR_EXCEPTION
299 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
300 "setParameters: Parameter \"partitioner: local parts\" "
301 "has the wrong type.");
305 if (List.isParameter (
"partitioner: overlap level")) {
306 if (List.isType<
int> (
"partitioner: overlap level")) {
307 OverlapLevel_ = List.get<
int> (
"partitioner: overlap level");
309 else if (! std::is_same<int, local_ordinal_type>::value &&
314 TEUCHOS_TEST_FOR_EXCEPTION
315 (
true, std::invalid_argument,
"Ifpack2::BlockRelaxation::"
316 "setParameters: Parameter \"partitioner: overlap level\" "
317 "has the wrong type.");
322 if ( ( List.isParameter(
"partitioner: global ID parts")) && (OverlapLevel_ < 1)) OverlapLevel_ = 1;
324 if (List.isParameter (
"partitioner: nonsymmetric overlap combine"))
325 nonsymCombine_ = List.get<
bool> (
"partitioner: nonsymmetric overlap combine");
327 if (List.isParameter (
"partitioner: combine mode"))
328 schwarzCombineMode_ = List.get<std::string> (
"partitioner: combine mode");
330 std::string defaultContainer =
"TriDi";
337 if (PrecType_ != Ifpack2::Details::JACOBI) {
341 NumLocalBlocks_ = A_->getLocalNumRows() / (-NumLocalBlocks_);
345 if(List.isParameter(
"block relaxation: decouple dofs"))
346 decouple_ = List.get<
bool>(
"block relaxation: decouple dofs");
351 TEUCHOS_TEST_FOR_EXCEPTION(
352 DoBackwardGS_, std::runtime_error,
353 "Ifpack2::BlockRelaxation:setParameters: Setting the \"relaxation: "
354 "backward mode\" parameter to true is not yet supported.");
356 if(List.isParameter(
"timer for apply"))
357 TimerForApply_ = List.get<
bool>(
"timer for apply");
364template<
class MatrixType,
class ContainerType>
365Teuchos::RCP<const Teuchos::Comm<int> >
368 TEUCHOS_TEST_FOR_EXCEPTION
369 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::getComm: "
370 "The matrix is null. You must call setMatrix() with a nonnull matrix "
371 "before you may call this method.");
372 return A_->getComm ();
375template<
class MatrixType,
class ContainerType>
376Teuchos::RCP<
const Tpetra::RowMatrix<
typename MatrixType::scalar_type,
377 typename MatrixType::local_ordinal_type,
378 typename MatrixType::global_ordinal_type,
379 typename MatrixType::node_type> >
384template<
class MatrixType,
class ContainerType>
385Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
386 typename MatrixType::global_ordinal_type,
387 typename MatrixType::node_type> >
391 TEUCHOS_TEST_FOR_EXCEPTION
392 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
393 "getDomainMap: The matrix is null. You must call setMatrix() with a "
394 "nonnull matrix before you may call this method.");
395 return A_->getDomainMap ();
398template<
class MatrixType,
class ContainerType>
399Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
400 typename MatrixType::global_ordinal_type,
401 typename MatrixType::node_type> >
405 TEUCHOS_TEST_FOR_EXCEPTION
406 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
407 "getRangeMap: The matrix is null. You must call setMatrix() with a "
408 "nonnull matrix before you may call this method.");
409 return A_->getRangeMap ();
412template<
class MatrixType,
class ContainerType>
419template<
class MatrixType,
class ContainerType>
423 return NumInitialize_;
426template<
class MatrixType,
class ContainerType>
434template<
class MatrixType,
class ContainerType>
442template<
class MatrixType,
class ContainerType>
447 return InitializeTime_;
450template<
class MatrixType,
class ContainerType>
458template<
class MatrixType,
class ContainerType>
467template<
class MatrixType,
class ContainerType>
469 TEUCHOS_TEST_FOR_EXCEPTION(
470 A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::getNodeSmootherComplexity: "
471 "The input matrix A is null. Please call setMatrix() with a nonnull "
472 "input matrix, then call compute(), before calling this method.");
475 size_t block_nnz = 0;
477 block_nnz += Partitioner_->numRowsInPart(i) *Partitioner_->numRowsInPart(i);
479 return block_nnz + A_->getLocalNumEntries();
482template<
class MatrixType,
class ContainerType>
485apply (
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
486 typename MatrixType::local_ordinal_type,
487 typename MatrixType::global_ordinal_type,
488 typename MatrixType::node_type>& X,
489 Tpetra::MultiVector<
typename MatrixType::scalar_type,
490 typename MatrixType::local_ordinal_type,
491 typename MatrixType::global_ordinal_type,
492 typename MatrixType::node_type>& Y,
493 Teuchos::ETransp mode,
497 TEUCHOS_TEST_FOR_EXCEPTION
498 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::apply: "
499 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
500 "then call initialize() and compute() (in that order), before you may "
501 "call this method.");
502 TEUCHOS_TEST_FOR_EXCEPTION(
503 !
isComputed (), std::runtime_error,
"Ifpack2::BlockRelaxation::apply: "
504 "isComputed() must be true prior to calling apply.");
505 TEUCHOS_TEST_FOR_EXCEPTION(
506 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
507 "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = "
508 << X.getNumVectors () <<
" != Y.getNumVectors() = "
509 << Y.getNumVectors () <<
".");
510 TEUCHOS_TEST_FOR_EXCEPTION(
511 mode != Teuchos::NO_TRANS, std::logic_error,
"Ifpack2::BlockRelaxation::"
512 "apply: This method currently only implements the case mode == Teuchos::"
513 "NO_TRANS. Transposed modes are not currently supported.");
514 TEUCHOS_TEST_FOR_EXCEPTION(
515 alpha != Teuchos::ScalarTraits<scalar_type>::one (), std::logic_error,
516 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
517 "the case alpha == 1. You specified alpha = " << alpha <<
".");
518 TEUCHOS_TEST_FOR_EXCEPTION(
519 beta != Teuchos::ScalarTraits<scalar_type>::zero (), std::logic_error,
520 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
521 "the case beta == 0. You specified beta = " << beta <<
".");
523 const std::string timerName (
"Ifpack2::BlockRelaxation::apply");
524 Teuchos::RCP<Teuchos::Time> timer;
525 if (TimerForApply_) {
526 timer = Teuchos::TimeMonitor::lookupCounter (timerName);
527 if (timer.is_null ()) {
528 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
532 Teuchos::Time time = Teuchos::Time(timerName);
533 double startTime = time.wallTime();
536 Teuchos::RCP<Teuchos::TimeMonitor> timeMon;
538 timeMon = Teuchos::rcp(
new Teuchos::TimeMonitor(*timer));
542 Teuchos::RCP<const MV> X_copy;
545 X_copy = rcp (
new MV (X, Teuchos::Copy));
547 X_copy = rcpFromRef (X);
551 if (ZeroStartingSolution_) {
552 Y.putScalar (STS::zero ());
557 case Ifpack2::Details::JACOBI:
558 ApplyInverseJacobi(*X_copy,Y);
560 case Ifpack2::Details::GS:
561 ApplyInverseGS(*X_copy,Y);
563 case Ifpack2::Details::SGS:
564 ApplyInverseSGS(*X_copy,Y);
566 case Ifpack2::Details::MTSPLITJACOBI:
568 Container_->applyInverseJacobi(*X_copy, Y, DampingFactor_, ZeroStartingSolution_, NumSweeps_);
571 TEUCHOS_TEST_FOR_EXCEPTION
572 (
true, std::logic_error,
"Ifpack2::BlockRelaxation::apply: Invalid "
573 "PrecType_ enum value " << PrecType_ <<
". Valid values are Ifpack2::"
574 "Details::JACOBI = " << Ifpack2::Details::JACOBI <<
", Ifpack2::Details"
575 "::GS = " << Ifpack2::Details::GS <<
", and Ifpack2::Details::SGS = "
576 << Ifpack2::Details::SGS <<
". Please report this bug to the Ifpack2 "
581 ApplyTime_ += (time.wallTime() - startTime);
585template<
class MatrixType,
class ContainerType>
588applyMat (
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
589 typename MatrixType::local_ordinal_type,
590 typename MatrixType::global_ordinal_type,
591 typename MatrixType::node_type>& X,
592 Tpetra::MultiVector<
typename MatrixType::scalar_type,
593 typename MatrixType::local_ordinal_type,
594 typename MatrixType::global_ordinal_type,
595 typename MatrixType::node_type>& Y,
596 Teuchos::ETransp mode)
const
598 A_->apply (X, Y, mode);
601template<
class MatrixType,
class ContainerType>
607 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
610 TEUCHOS_TEST_FOR_EXCEPTION
611 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::initialize: "
612 "The matrix is null. You must call setMatrix() with a nonnull matrix "
613 "before you may call this method.");
615 Teuchos::RCP<Teuchos::Time> timer =
616 Teuchos::TimeMonitor::getNewCounter (
"Ifpack2::BlockRelaxation::initialize");
617 double startTime = timer->wallTime();
620 Teuchos::TimeMonitor timeMon (*timer);
621 IsInitialized_ =
false;
624 Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
625 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A_);
626 hasBlockCrsMatrix_ = !A_bcrs.is_null();
627 if (A_bcrs.is_null ()) {
628 hasBlockCrsMatrix_ =
false;
631 hasBlockCrsMatrix_ =
true;
634 NumLocalRows_ = A_->getLocalNumRows ();
635 NumGlobalRows_ = A_->getGlobalNumRows ();
636 NumGlobalNonzeros_ = A_->getGlobalNumEntries ();
641 Partitioner_ = Teuchos::null;
643 if (PartitionerType_ ==
"linear") {
646 }
else if (PartitionerType_ ==
"line") {
649 }
else if (PartitionerType_ ==
"user") {
655 TEUCHOS_TEST_FOR_EXCEPTION
656 (
true, std::logic_error,
"Ifpack2::BlockRelaxation::initialize: Unknown "
657 "partitioner type " << PartitionerType_ <<
". Valid values are "
658 "\"linear\", \"line\", and \"user\".");
662 Partitioner_->setParameters (List_);
663 Partitioner_->compute ();
666 NumLocalBlocks_ = Partitioner_->numLocalParts ();
671 if (A_->getComm()->getSize() != 1) {
679 TEUCHOS_TEST_FOR_EXCEPTION
680 (NumSweeps_ < 0, std::logic_error,
"Ifpack2::BlockRelaxation::initialize: "
681 "NumSweeps_ = " << NumSweeps_ <<
" < 0.");
684 ExtractSubmatricesStructure ();
688 if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) {
689 TEUCHOS_TEST_FOR_EXCEPTION
690 (hasBlockCrsMatrix_, std::runtime_error,
691 "Ifpack2::BlockRelaxation::initialize: "
692 "We do not support overlapped Jacobi yet for Tpetra::BlockCrsMatrix. Sorry!");
695 W_ = rcp (
new vector_type (A_->getRowMap ()));
696 W_->putScalar (STS::zero ());
698 Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0);
701 for (
size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) {
703 w_ptr[LID] += STS::one();
711 if (schwarzCombineMode_ ==
"ADD") {
712 typedef Tpetra::MultiVector< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type,typename MatrixType::node_type> scMV;
713 Teuchos::RCP<const import_type> theImport = A_->getGraph()->getImporter();
714 if (!theImport.is_null()) {
715 scMV nonOverLapW(theImport->getSourceMap(), 1,
false);
716 Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0);
717 Teuchos::ArrayRCP<scalar_type> nonOverLapWArray = nonOverLapW.getDataNonConst(0);
718 nonOverLapW.putScalar(STS::zero ());
719 for (
int ii = 0; ii < (int) theImport->getSourceMap()->getLocalNumElements(); ii++) nonOverLapWArray[ii] = w_ptr[ii];
720 nonOverLapWArray = Teuchos::null;
721 w_ptr = Teuchos::null;
722 nonOverLapW.doExport (*W_, *theImport, Tpetra::ADD);
723 W_->doImport( nonOverLapW, *theImport, Tpetra::INSERT);
727 W_->reciprocal (*W_);
731 InitializeTime_ += (timer->wallTime() - startTime);
733 IsInitialized_ =
true;
737template<
class MatrixType,
class ContainerType>
744 TEUCHOS_TEST_FOR_EXCEPTION
745 (A_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::compute: "
746 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
747 "then call initialize(), before you may call this method.");
753 Teuchos::RCP<Teuchos::Time> timer =
754 Teuchos::TimeMonitor::getNewCounter (
"Ifpack2::BlockRelaxation::compute");
756 double startTime = timer->wallTime();
759 Teuchos::TimeMonitor timeMon (*timer);
764 Container_->compute();
767 ComputeTime_ += (timer->wallTime() - startTime);
772template<
class MatrixType,
class ContainerType>
777 TEUCHOS_TEST_FOR_EXCEPTION
778 (Partitioner_.is_null (), std::runtime_error,
"Ifpack2::BlockRelaxation::"
779 "ExtractSubmatricesStructure: Partitioner object is null.");
781 std::string containerType = ContainerType::getName ();
782 if (containerType ==
"Generic") {
786 containerType = containerType_;
790 bool pointIndexed = decouple_ && hasBlockCrsMatrix_;
791 Teuchos::Array<Teuchos::Array<local_ordinal_type> > blockRows;
795 auto A_bcrs = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(A_);
796 local_ordinal_type dofs = hasBlockCrsMatrix_ ?
797 A_bcrs->getBlockSize() : List_.get<
int>(
"partitioner: PDE equations");
798 blockRows.resize(NumLocalBlocks_ * dofs);
799 if(hasBlockCrsMatrix_)
801 for(local_ordinal_type i = 0; i < NumLocalBlocks_; i++)
803 size_t blockSize = Partitioner_->numRowsInPart(i);
806 for(local_ordinal_type j = 0; j < dofs; j++)
808 local_ordinal_type blockIndex = i * dofs + j;
809 blockRows[blockIndex].resize(blockSize);
810 for(
size_t k = 0; k < blockSize; k++)
814 blockRows[blockIndex][k] = (*Partitioner_)(i, k) * dofs + j;
826 TEUCHOS_TEST_FOR_EXCEPTION(Partitioner_->numRowsInPart(i) % dofs != 0, std::logic_error,
827 "Expected size of all blocks (partitions) to be divisible by MueLu dofs/node.");
828 size_t blockSize = Partitioner_->numRowsInPart(i) / dofs;
835 blockRows[blockIndex].resize(blockSize);
836 for(
size_t k = 0; k < blockSize; k++)
838 blockRows[blockIndex][k] = (*Partitioner_)(i, k * dofs + j);
847 blockRows.resize(NumLocalBlocks_);
850 const size_t numRows = Partitioner_->numRowsInPart (i);
851 blockRows[i].resize(numRows);
853 for (
size_t j = 0; j < numRows; ++j)
855 blockRows[i][j] = (*Partitioner_) (i,j);
861 Container_->setParameters(List_);
862 Container_->initialize();
865template<
class MatrixType,
class ContainerType>
870 const size_t NumVectors = X.getNumVectors ();
872 MV AY (Y.getMap (), NumVectors);
875 int starting_iteration = 0;
876 if (OverlapLevel_ > 0)
879 auto WView = W_->getLocalViewHost (Tpetra::Access::ReadOnly);
880 if(ZeroStartingSolution_) {
881 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
882 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
883 Container_->DoOverlappingJacobi(XView, YView, WView, DampingFactor_, nonsymCombine_);
884 starting_iteration = 1;
886 const scalar_type ONE = STS::one();
887 for(
int j = starting_iteration; j < NumSweeps_; j++)
890 AY.update (ONE, X, -ONE);
892 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
893 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
894 Container_->DoOverlappingJacobi (AYView, YView, WView, DampingFactor_, nonsymCombine_);
901 if(ZeroStartingSolution_)
903 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
904 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
905 Container_->DoJacobi (XView, YView, DampingFactor_);
906 starting_iteration = 1;
908 const scalar_type ONE = STS::one();
909 for(
int j = starting_iteration; j < NumSweeps_; j++)
912 AY.update (ONE, X, -ONE);
914 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
915 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
916 Container_->DoJacobi (AYView, YView, DampingFactor_);
922template<
class MatrixType,
class ContainerType>
929 size_t numVecs = X.getNumVectors();
931 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
932 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
935 bool deleteY2 =
false;
938 Y2 = ptr(
new MV(Importer_->getTargetMap(), numVecs));
945 for(
int j = 0; j < NumSweeps_; ++j)
948 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
949 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
950 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
955 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
956 for(
int j = 0; j < NumSweeps_; ++j)
958 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
965template<
class MatrixType,
class ContainerType>
973 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
974 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
977 bool deleteY2 =
false;
980 Y2 = ptr(
new MV(Importer_->getTargetMap(), X.getNumVectors()));
987 for(
int j = 0; j < NumSweeps_; ++j)
990 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
991 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
992 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
997 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
998 for(
int j = 0; j < NumSweeps_; ++j)
1000 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
1007template<
class MatrixType,
class ContainerType>
1008void BlockRelaxation<MatrixType,ContainerType>::computeImporter ()
const
1013 using Teuchos::ArrayView;
1014 using Teuchos::Array;
1015 using Teuchos::Comm;
1016 using Teuchos::rcp_dynamic_cast;
1019 if(hasBlockCrsMatrix_)
1021 const RCP<const block_crs_matrix_type> bcrs =
1022 rcp_dynamic_cast<const block_crs_matrix_type>(A_);
1023 int bs_ = bcrs->getBlockSize();
1024 RCP<const map_type> oldDomainMap = A_->getDomainMap();
1025 RCP<const map_type> oldColMap = A_->getColMap();
1028 global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
1030 RCP<const Comm<int>> comm = oldColMap->getComm();
1031 ArrayView<const global_ordinal_type> oldColElements = oldColMap->getLocalElementList();
1032 Array<global_ordinal_type> newColElements(bs_ * oldColElements.size());
1033 for(
int i = 0; i < oldColElements.size(); i++)
1035 for(
int j = 0; j < bs_; j++)
1036 newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
1038 RCP<map_type> colMap(
new map_type(numGlobalElements, newColElements, indexBase, comm));
1040 Importer_ = rcp(
new import_type(oldDomainMap, colMap));
1042 else if(!A_.is_null())
1044 Importer_ = A_->getGraph()->getImporter();
1045 if(Importer_.is_null())
1046 Importer_ = rcp(
new import_type(A_->getDomainMap(), A_->getColMap()));
1052template<
class MatrixType,
class ContainerType>
1057 std::ostringstream out;
1062 out <<
"\"Ifpack2::BlockRelaxation\": {";
1063 if (this->getObjectLabel () !=
"") {
1064 out <<
"Label: \"" << this->getObjectLabel () <<
"\", ";
1068 if (A_.is_null ()) {
1069 out <<
"Matrix: null, ";
1074 out <<
"Global matrix dimensions: ["
1075 << A_->getGlobalNumRows () <<
", " << A_->getGlobalNumCols () <<
"], ";
1080 out <<
"\"relaxation: type\": ";
1081 if (PrecType_ == Ifpack2::Details::JACOBI) {
1082 out <<
"Block Jacobi";
1083 }
else if (PrecType_ == Ifpack2::Details::GS) {
1084 out <<
"Block Gauss-Seidel";
1085 }
else if (PrecType_ == Ifpack2::Details::SGS) {
1086 out <<
"Block Symmetric Gauss-Seidel";
1087 }
else if (PrecType_ == Ifpack2::Details::MTSPLITJACOBI) {
1088 out <<
"MT Split Jacobi";
1094 if(hasBlockCrsMatrix_)
1098 int approx_rows_per_part = A_->getLocalNumRows()/Partitioner_->numLocalParts();
1099 out <<
", blocksize: "<<approx_rows_per_part;
1101 out <<
", overlap: " << OverlapLevel_;
1103 out <<
", " <<
"sweeps: " << NumSweeps_ <<
", "
1104 <<
"damping factor: " << DampingFactor_ <<
", ";
1106 std::string containerType = ContainerType::getName();
1107 out <<
"block container: " << ((containerType ==
"Generic") ? containerType_ : containerType);
1108 if(List_.isParameter(
"partitioner: PDE equations"))
1109 out <<
", dofs/node: "<<List_.get<
int>(
"partitioner: PDE equations");
1116template<
class MatrixType,
class ContainerType>
1119describe (Teuchos::FancyOStream& out,
1120 const Teuchos::EVerbosityLevel verbLevel)
const
1124 using Teuchos::TypeNameTraits;
1125 using Teuchos::VERB_DEFAULT;
1126 using Teuchos::VERB_NONE;
1127 using Teuchos::VERB_LOW;
1128 using Teuchos::VERB_MEDIUM;
1129 using Teuchos::VERB_HIGH;
1130 using Teuchos::VERB_EXTREME;
1132 Teuchos::EVerbosityLevel vl = verbLevel;
1133 if (vl == VERB_DEFAULT) vl = VERB_LOW;
1134 const int myImageID = A_->getComm()->getRank();
1137 Teuchos::OSTab tab (out);
1144 if (vl != VERB_NONE && myImageID == 0) {
1145 out <<
"Ifpack2::BlockRelaxation<"
1146 << TypeNameTraits<MatrixType>::name () <<
", "
1147 << TypeNameTraits<ContainerType>::name () <<
" >:";
1148 Teuchos::OSTab tab1 (out);
1150 if (this->getObjectLabel () !=
"") {
1151 out <<
"label: \"" << this->getObjectLabel () <<
"\"" << endl;
1153 out <<
"initialized: " << (
isInitialized () ?
"true" :
"false") << endl
1154 <<
"computed: " << (
isComputed () ?
"true" :
"false") << endl;
1156 std::string precType;
1157 if (PrecType_ == Ifpack2::Details::JACOBI) {
1158 precType =
"Block Jacobi";
1159 }
else if (PrecType_ == Ifpack2::Details::GS) {
1160 precType =
"Block Gauss-Seidel";
1161 }
else if (PrecType_ == Ifpack2::Details::SGS) {
1162 precType =
"Block Symmetric Gauss-Seidel";
1164 out <<
"type: " << precType << endl
1165 <<
"global number of rows: " << A_->getGlobalNumRows () << endl
1166 <<
"global number of columns: " << A_->getGlobalNumCols () << endl
1167 <<
"number of sweeps: " << NumSweeps_ << endl
1168 <<
"damping factor: " << DampingFactor_ << endl
1169 <<
"nonsymmetric overlap combine" << nonsymCombine_ << endl
1170 <<
"backwards mode: "
1171 << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ?
"true" :
"false")
1173 <<
"zero starting solution: "
1174 << (ZeroStartingSolution_ ?
"true" :
"false") << endl;
1192#ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1194#define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \
1196 class Ifpack2::BlockRelaxation< \
1197 Tpetra::RowMatrix<S, LO, GO, N>, \
1198 Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N> > >;
Ifpack2::BlockRelaxation class declaration.
Declaration of a user-defined partitioner in which the user can define a partition of the graph....
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition Ifpack2_BlockRelaxation_decl.hpp:91
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner's constructor.
Definition Ifpack2_BlockRelaxation_def.hpp:380
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:389
int getNumApply() const
Returns the number of calls to apply().
Definition Ifpack2_BlockRelaxation_def.hpp:437
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition Ifpack2_BlockRelaxation_def.hpp:127
BlockRelaxation(const Teuchos::RCP< const row_matrix_type > &Matrix)
Constructor.
Definition Ifpack2_BlockRelaxation_def.hpp:86
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition Ifpack2_BlockRelaxation_def.hpp:740
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_BlockRelaxation_def.hpp:468
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:1119
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:100
virtual ~BlockRelaxation()
Destructor.
Definition Ifpack2_BlockRelaxation_def.hpp:121
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition Ifpack2_BlockRelaxation_decl.hpp:231
int getNumCompute() const
Returns the number of calls to compute().
Definition Ifpack2_BlockRelaxation_def.hpp:429
double getApplyTime() const
Returns the time spent in apply().
Definition Ifpack2_BlockRelaxation_def.hpp:461
int getNumInitialize() const
Returns the number of calls to initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:422
void initialize()
Initialize.
Definition Ifpack2_BlockRelaxation_def.hpp:604
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:403
bool isComputed() const
Return true if compute() has been called.
Definition Ifpack2_BlockRelaxation_decl.hpp:239
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:588
double getComputeTime() const
Returns the time spent in compute().
Definition Ifpack2_BlockRelaxation_def.hpp:453
void setParameters(const Teuchos::ParameterList ¶ms)
Sets all the parameters for the preconditioner.
Definition Ifpack2_BlockRelaxation_def.hpp:181
double getInitializeTime() const
Returns the time spent in initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:445
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:485
std::string description() const
A one-line description of this object.
Definition Ifpack2_BlockRelaxation_def.hpp:1055
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition Ifpack2_BlockRelaxation_def.hpp:366
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_BlockRelaxation_def.hpp:58
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:97
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
Ifpack2::LinePartitioner: A class to define partitions into a set of lines.
Definition Ifpack2_LinePartitioner_decl.hpp:78
A class to define linear partitions.
Definition Ifpack2_LinearPartitioner_decl.hpp:60
Ifpack2 implementation details.
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
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
static Teuchos::RCP< BaseContainer > build(std::string containerType, const Teuchos::RCP< const MatrixType > &A, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, const Teuchos::RCP< const import_type > importer, bool pointIndexed)
Build a specialization of Ifpack2::Container given a key that has been registered.
Definition Ifpack2_ContainerFactory_def.hpp:89