10 #ifndef IFPACK2_BLOCKRELAXATION_DEF_HPP
11 #define IFPACK2_BLOCKRELAXATION_DEF_HPP
14 #include "Ifpack2_LinearPartitioner.hpp"
15 #include "Ifpack2_LinePartitioner.hpp"
16 #include "Ifpack2_Zoltan2Partitioner_decl.hpp"
17 #include "Ifpack2_Zoltan2Partitioner_def.hpp"
19 #include "Ifpack2_Details_UserPartitioner_def.hpp"
20 #include "Ifpack2_LocalFilter.hpp"
21 #include "Ifpack2_Parameters.hpp"
23 #include "Tpetra_BlockCrsMatrix_Helpers_decl.hpp"
24 #include "Tpetra_Import_Util.hpp"
25 #include "Ifpack2_BlockHelper_Timers.hpp"
29 template <
class MatrixType,
class ContainerType>
33 IsInitialized_ =
false;
35 Partitioner_ = Teuchos::null;
39 IsParallel_ = (A->getRowMap()->getComm()->getSize() != 1);
41 Teuchos::rcp_dynamic_cast<
const block_crs_matrix_type>(A);
42 hasBlockCrsMatrix_ = !A_bcrs.
is_null();
44 if (!Container_.is_null()) {
47 Container_->clearBlocks();
56 template <
class MatrixType,
class ContainerType>
59 : Container_(Teuchos::null)
60 , Partitioner_(Teuchos::null)
61 , PartitionerType_(
"linear")
64 , containerType_(
"TriDi")
65 , PrecType_(Ifpack2::Details::JACOBI)
66 , ZeroStartingSolution_(true)
67 , hasBlockCrsMatrix_(false)
68 , DoBackwardGS_(false)
71 , schwarzCombineMode_(
"ZERO")
72 , DampingFactor_(
STS::one())
73 , IsInitialized_(false)
77 , TimerForApply_(true)
79 , InitializeTime_(0.0)
83 , NumGlobalNonzeros_(0)
85 , Importer_(Teuchos::null) {
89 template <
class MatrixType,
class ContainerType>
93 template <
class MatrixType,
class ContainerType>
99 validParams->
set(
"relaxation: container",
"TriDi");
100 validParams->
set(
"relaxation: backward mode",
false);
101 validParams->
set(
"relaxation: type",
"Jacobi");
102 validParams->
set(
"relaxation: sweeps", 1);
103 validParams->
set(
"relaxation: damping factor", STS::one());
104 validParams->
set(
"relaxation: zero starting solution",
true);
105 validParams->
set(
"block relaxation: decouple dofs",
false);
106 validParams->
set(
"schwarz: compute condest",
false);
107 validParams->
set(
"schwarz: combine mode",
"ZERO");
108 validParams->
set(
"schwarz: use reordering",
true);
109 validParams->
set(
"schwarz: filter singletons",
false);
110 validParams->
set(
"schwarz: overlap level", 0);
111 validParams->
set(
"partitioner: type",
"greedy");
112 validParams->
set(
"zoltan2: algorithm",
"phg");
113 validParams->
set(
"partitioner: local parts", 1);
114 validParams->
set(
"partitioner: overlap", 0);
115 validParams->
set(
"partitioner: combine mode",
"ZERO");
117 validParams->
set(
"partitioner: parts", tmp0);
119 validParams->
set(
"partitioner: global ID parts", tmp1);
120 validParams->
set(
"partitioner: nonsymmetric overlap combine",
false);
121 validParams->
set(
"partitioner: maintain sparsity",
false);
122 validParams->
set(
"fact: ilut level-of-fill", 1.0);
123 validParams->
set(
"fact: absolute threshold", 0.0);
124 validParams->
set(
"fact: relative threshold", 1.0);
125 validParams->
set(
"fact: relax value", 0.0);
128 validParams->
set(
"Amesos2", dummyList);
130 validParams->
set(
"Amesos2 solver name",
"KLU2");
133 validParams->
set(
"partitioner: map", tmp);
135 validParams->
set(
"partitioner: line detection threshold", 0.0);
136 validParams->
set(
"partitioner: PDE equations", 1);
138 typename MatrixType::local_ordinal_type,
139 typename MatrixType::global_ordinal_type,
140 typename MatrixType::node_type>>
142 validParams->
set(
"partitioner: coordinates", dummy);
143 validParams->
set(
"timer for apply",
true);
144 validParams->
set(
"partitioner: subparts per part", 1);
145 validParams->
set(
"partitioner: block size", -1);
146 validParams->
set(
"partitioner: print level",
false);
147 validParams->
set(
"partitioner: explicit convert to BlockCrs",
false);
148 validParams->
set(
"partitioner: checkBlockConsistency",
true);
149 validParams->
set(
"partitioner: use LIDs",
true);
154 template <
class MatrixType,
class ContainerType>
160 this->setParametersImpl(const_cast<Teuchos::ParameterList&>(pl));
163 template <
class MatrixType,
class ContainerType>
166 if (List.
isType<
double>(
"relaxation: damping factor")) {
169 scalar_type df = List.
get<
double>(
"relaxation: damping factor");
170 List.
remove(
"relaxation: damping factor");
171 List.
set(
"relaxation: damping factor", df);
185 containerType_ = List.
get<std::string>(
"relaxation: container");
187 if (containerType_ ==
"Tridiagonal") {
188 containerType_ =
"TriDi";
190 if (containerType_ ==
"Block Tridiagonal") {
191 containerType_ =
"BlockTriDi";
196 const std::string relaxType = List.
get<std::string>(
"relaxation: type");
198 if (relaxType ==
"Jacobi") {
199 PrecType_ = Ifpack2::Details::JACOBI;
200 }
else if (relaxType ==
"MT Split Jacobi") {
201 PrecType_ = Ifpack2::Details::MTSPLITJACOBI;
202 }
else if (relaxType ==
"Gauss-Seidel") {
203 PrecType_ = Ifpack2::Details::GS;
204 }
else if (relaxType ==
"Symmetric Gauss-Seidel") {
205 PrecType_ = Ifpack2::Details::SGS;
208 "Ifpack2::BlockRelaxation::"
209 "setParameters: Invalid parameter value \""
211 <<
"\" for parameter \"relaxation: type\".");
216 NumSweeps_ = List.
get<
int>(
"relaxation: sweeps");
222 if (List.
isParameter(
"relaxation: damping factor")) {
223 if (List.
isType<
double>(
"relaxation: damping factor")) {
224 const double dampingFactor =
225 List.
get<
double>(
"relaxation: damping factor");
226 DampingFactor_ =
static_cast<scalar_type
>(dampingFactor);
227 }
else if (List.
isType<scalar_type>(
"relaxation: damping factor")) {
228 DampingFactor_ = List.
get<scalar_type>(
"relaxation: damping factor");
229 }
else if (List.
isType<magnitude_type>(
"relaxation: damping factor")) {
230 const magnitude_type dampingFactor =
231 List.
get<magnitude_type>(
"relaxation: damping factor");
232 DampingFactor_ =
static_cast<scalar_type
>(dampingFactor);
235 "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")) {
256 if (List.
isType<local_ordinal_type>(
"partitioner: local parts")) {
257 NumLocalBlocks_ = List.
get<local_ordinal_type>(
"partitioner: local parts");
258 }
else if (!std::is_same<int, local_ordinal_type>::value &&
259 List.
isType<
int>(
"partitioner: local parts")) {
260 NumLocalBlocks_ = List.
get<
int>(
"partitioner: local parts");
263 "Ifpack2::BlockRelaxation::"
264 "setParameters: Parameter \"partitioner: local parts\" "
265 "has the wrong type.");
269 if (List.
isParameter(
"partitioner: overlap level")) {
270 if (List.
isType<
int>(
"partitioner: overlap level")) {
271 OverlapLevel_ = List.
get<
int>(
"partitioner: overlap level");
272 }
else if (!std::is_same<int, local_ordinal_type>::value &&
273 List.
isType<local_ordinal_type>(
"partitioner: overlap level")) {
274 OverlapLevel_ = List.
get<local_ordinal_type>(
"partitioner: overlap level");
277 "Ifpack2::BlockRelaxation::"
278 "setParameters: Parameter \"partitioner: overlap level\" "
279 "has the wrong type.");
284 if ((List.
isParameter(
"partitioner: global ID parts")) && (OverlapLevel_ < 1)) OverlapLevel_ = 1;
286 if (List.
isParameter(
"partitioner: nonsymmetric overlap combine"))
287 nonsymCombine_ = List.
get<
bool>(
"partitioner: nonsymmetric overlap combine");
290 schwarzCombineMode_ = List.
get<std::string>(
"partitioner: combine mode");
292 std::string defaultContainer =
"TriDi";
293 if (std::is_same<ContainerType, Container<MatrixType>>::value) {
298 if (PrecType_ != Ifpack2::Details::JACOBI) {
301 if (NumLocalBlocks_ < static_cast<local_ordinal_type>(0)) {
302 NumLocalBlocks_ = A_->getLocalNumRows() / (-NumLocalBlocks_);
306 if (List.
isParameter(
"block relaxation: decouple dofs"))
307 decouple_ = List.
get<
bool>(
"block relaxation: decouple dofs");
313 DoBackwardGS_, std::runtime_error,
314 "Ifpack2::BlockRelaxation:setParameters: Setting the \"relaxation: "
315 "backward mode\" parameter to true is not yet supported.");
318 TimerForApply_ = List.
get<
bool>(
"timer for apply");
325 template <
class MatrixType,
class ContainerType>
329 "Ifpack2::BlockRelaxation::getComm: "
330 "The matrix is null. You must call setMatrix() with a nonnull matrix "
331 "before you may call this method.");
332 return A_->getComm();
335 template <
class MatrixType,
class ContainerType>
336 Teuchos::RCP<
const Tpetra::RowMatrix<
typename MatrixType::scalar_type,
337 typename MatrixType::local_ordinal_type,
338 typename MatrixType::global_ordinal_type,
339 typename MatrixType::node_type>>
344 template <
class MatrixType,
class ContainerType>
345 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
346 typename MatrixType::global_ordinal_type,
347 typename MatrixType::node_type>>
351 "Ifpack2::BlockRelaxation::"
352 "getDomainMap: The matrix is null. You must call setMatrix() with a "
353 "nonnull matrix before you may call this method.");
354 return A_->getDomainMap();
357 template <
class MatrixType,
class ContainerType>
358 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
359 typename MatrixType::global_ordinal_type,
360 typename MatrixType::node_type>>
364 "Ifpack2::BlockRelaxation::"
365 "getRangeMap: The matrix is null. You must call setMatrix() with a "
366 "nonnull matrix before you may call this method.");
367 return A_->getRangeMap();
370 template <
class MatrixType,
class ContainerType>
376 template <
class MatrixType,
class ContainerType>
379 return NumInitialize_;
382 template <
class MatrixType,
class ContainerType>
388 template <
class MatrixType,
class ContainerType>
394 template <
class MatrixType,
class ContainerType>
398 return InitializeTime_;
401 template <
class MatrixType,
class ContainerType>
408 template <
class MatrixType,
class ContainerType>
415 template <
class MatrixType,
class ContainerType>
418 A_.
is_null(), std::runtime_error,
419 "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_->getLocalNumEntries();
431 template <
class MatrixType,
class ContainerType>
433 apply(
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
434 typename MatrixType::local_ordinal_type,
435 typename MatrixType::global_ordinal_type,
436 typename MatrixType::node_type>& X,
437 Tpetra::MultiVector<
typename MatrixType::scalar_type,
438 typename MatrixType::local_ordinal_type,
439 typename MatrixType::global_ordinal_type,
440 typename MatrixType::node_type>& Y,
445 "Ifpack2::BlockRelaxation::apply: "
446 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
447 "then call initialize() and compute() (in that order), before you may "
448 "call this method.");
450 !isComputed(), std::runtime_error,
451 "Ifpack2::BlockRelaxation::apply: "
452 "isComputed() must be true prior to calling apply.");
453 TEUCHOS_TEST_FOR_EXCEPTION(
454 X.getNumVectors() != Y.getNumVectors(), std::invalid_argument,
455 "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = "
456 << X.getNumVectors() <<
" != Y.getNumVectors() = "
457 << Y.getNumVectors() <<
".");
458 TEUCHOS_TEST_FOR_EXCEPTION(
460 "Ifpack2::BlockRelaxation::"
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 = "
468 TEUCHOS_TEST_FOR_EXCEPTION(
470 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
471 "the case beta == 0. You specified beta = "
474 const std::string timerName(
"Ifpack2::BlockRelaxation::apply");
476 if (TimerForApply_) {
498 X_copy = rcpFromRef(X);
502 if (ZeroStartingSolution_) {
503 Y.putScalar(STS::zero());
508 case Ifpack2::Details::JACOBI:
509 ApplyInverseJacobi(*X_copy, Y);
511 case Ifpack2::Details::GS:
512 ApplyInverseGS(*X_copy, Y);
514 case Ifpack2::Details::SGS:
515 ApplyInverseSGS(*X_copy, Y);
517 case Ifpack2::Details::MTSPLITJACOBI:
519 Container_->applyInverseJacobi(*X_copy, Y, DampingFactor_, ZeroStartingSolution_, NumSweeps_);
522 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
523 "Ifpack2::BlockRelaxation::apply: Invalid "
524 "PrecType_ enum value "
525 << PrecType_ <<
". Valid values are Ifpack2::"
527 << Ifpack2::Details::JACOBI <<
", Ifpack2::Details"
529 << Ifpack2::Details::GS <<
", and Ifpack2::Details::SGS = "
530 << Ifpack2::Details::SGS <<
". Please report this bug to the Ifpack2 "
535 ApplyTime_ += (time.
wallTime() - startTime);
539 template <
class MatrixType,
class ContainerType>
541 applyMat(
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
542 typename MatrixType::local_ordinal_type,
543 typename MatrixType::global_ordinal_type,
544 typename MatrixType::node_type>& X,
545 Tpetra::MultiVector<
typename MatrixType::scalar_type,
546 typename MatrixType::local_ordinal_type,
547 typename MatrixType::global_ordinal_type,
548 typename MatrixType::node_type>& Y,
550 A_->apply(X, Y, mode);
553 template <
class MatrixType,
class ContainerType>
557 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
561 "Ifpack2::BlockRelaxation::initialize: "
562 "The matrix is null. You must call setMatrix() with a nonnull matrix "
563 "before you may call this method.");
567 double startTime = timer->wallTime();
571 IsInitialized_ =
false;
575 Teuchos::rcp_dynamic_cast<
const block_crs_matrix_type>(A_);
576 hasBlockCrsMatrix_ = !A_bcrs.
is_null();
580 if (!hasBlockCrsMatrix_ && List_.isParameter(
"relaxation: container") && List_.get<std::string>(
"relaxation: container") ==
"BlockTriDi") {
581 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::convertToBlockCrsMatrix", convertToBlockCrsMatrix);
582 int block_size = List_.get<
int>(
"partitioner: block size");
583 bool use_explicit_conversion = List_.isParameter(
"partitioner: explicit convert to BlockCrs") && List_.get<
bool>(
"partitioner: explicit convert to BlockCrs");
584 TEUCHOS_TEST_FOR_EXCEPT_MSG(use_explicit_conversion && block_size == -1,
"A pointwise matrix and block_size = -1 were given as inputs.");
585 bool use_LID = !List_.isParameter(
"partitioner: use LIDs") || List_.get<
bool>(
"partitioner: use LIDs");
586 bool check_block_consistency = !List_.isParameter(
"partitioner: checkBlockConsistency") || List_.get<
bool>(
"partitioner: checkBlockConsistency");
588 if ((use_LID || !use_explicit_conversion) && check_block_consistency) {
589 if (!A_->getGraph()->getImporter().
is_null()) {
590 TEUCHOS_TEST_FOR_EXCEPT_MSG(!Tpetra::Import_Util::checkBlockConsistency(*(A_->getGraph()->getColMap()), block_size),
591 "The pointwise graph of the input matrix A pointwise is not consistent with block_size.");
594 if (use_explicit_conversion) {
595 A_bcrs = Tpetra::convertToBlockCrsMatrix(*Teuchos::rcp_dynamic_cast<const crs_matrix_type>(A_), block_size, use_LID);
597 hasBlockCrsMatrix_ =
true;
598 graph = A_->getGraph();
600 graph = Tpetra::getBlockCrsGraph(*Teuchos::rcp_dynamic_cast<const crs_matrix_type>(A_), block_size,
true);
602 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
605 NumLocalRows_ = A_->getLocalNumRows();
606 NumGlobalRows_ = A_->getGlobalNumRows();
607 NumGlobalNonzeros_ = A_->getGlobalNumEntries();
612 Partitioner_ = Teuchos::null;
614 if (PartitionerType_ ==
"linear") {
615 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::linear", linear);
618 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
619 }
else if (PartitionerType_ ==
"line") {
620 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::line", line);
623 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
624 }
else if (PartitionerType_ ==
"user") {
625 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::user", user);
628 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
629 }
else if (PartitionerType_ ==
"zoltan2") {
630 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::zoltan2", zoltan2);
631 #if defined(HAVE_IFPACK2_ZOLTAN2)
632 if (graph->getComm()->getSize() == 1) {
635 rcp(
new Ifpack2::Zoltan2Partitioner<row_graph_type>(graph));
640 rcp(
new Ifpack2::Zoltan2Partitioner<row_graph_type>(A_local->getGraph()));
645 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
650 "Ifpack2::BlockRelaxation::initialize: Unknown "
652 << PartitionerType_ <<
". Valid values are "
653 "\"linear\", \"line\", and \"user\".");
658 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::Partitioner",
Partitioner);
659 Partitioner_->setParameters(List_);
660 Partitioner_->compute();
661 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
665 NumLocalBlocks_ = Partitioner_->numLocalParts();
670 if (A_->getComm()->getSize() != 1) {
679 "Ifpack2::BlockRelaxation::initialize: "
681 << NumSweeps_ <<
" < 0.");
685 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::ExtractSubmatricesStructure", ExtractSubmatricesStructure);
686 ExtractSubmatricesStructure();
687 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
692 if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) {
694 "Ifpack2::BlockRelaxation::initialize: "
695 "We do not support overlapped Jacobi yet for Tpetra::BlockCrsMatrix. Sorry!");
698 W_ =
rcp(
new vector_type(A_->getRowMap()));
699 W_->putScalar(STS::zero());
704 for (
size_t j = 0; j < Partitioner_->numRowsInPart(i); ++j) {
706 w_ptr[LID] += STS::one();
714 if (schwarzCombineMode_ ==
"ADD") {
715 IFPACK2_BLOCKHELPER_TIMER(
"Ifpack2::BlockRelaxation::initialize::ADD", ADD);
716 typedef Tpetra::MultiVector<typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> scMV;
719 scMV nonOverLapW(theImport->getSourceMap(), 1,
false);
722 nonOverLapW.putScalar(STS::zero());
723 for (
int ii = 0; ii < (int)theImport->getSourceMap()->getLocalNumElements(); ii++) nonOverLapWArray[ii] = w_ptr[ii];
724 nonOverLapWArray = Teuchos::null;
725 w_ptr = Teuchos::null;
726 nonOverLapW.doExport(*W_, *theImport, Tpetra::ADD);
727 W_->doImport(nonOverLapW, *theImport, Tpetra::INSERT);
729 IFPACK2_BLOCKHELPER_TIMER_DEFAULT_FENCE();
735 InitializeTime_ += (timer->wallTime() - startTime);
737 IsInitialized_ =
true;
740 template <
class MatrixType,
class ContainerType>
746 "Ifpack2::BlockRelaxation::compute: "
747 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
748 "then call initialize(), before you may call this method.");
750 if (!isInitialized()) {
757 double startTime = timer->
wallTime();
765 Container_->compute();
768 ComputeTime_ += (timer->
wallTime() - startTime);
773 template <
class MatrixType,
class ContainerType>
777 "Ifpack2::BlockRelaxation::"
778 "ExtractSubmatricesStructure: Partitioner object is null.");
780 std::string containerType = ContainerType::getName();
781 if (containerType ==
"Generic") {
785 containerType = containerType_;
789 bool pointIndexed = decouple_ && hasBlockCrsMatrix_;
793 auto A_bcrs = Teuchos::rcp_dynamic_cast<
const block_crs_matrix_type>(A_);
794 local_ordinal_type dofs = hasBlockCrsMatrix_ ? A_bcrs->getBlockSize() : List_.get<
int>(
"partitioner: PDE equations");
795 blockRows.
resize(NumLocalBlocks_ * dofs);
796 if (hasBlockCrsMatrix_) {
797 for (local_ordinal_type i = 0; i < NumLocalBlocks_; i++) {
798 size_t blockSize = Partitioner_->numRowsInPart(i);
801 for (local_ordinal_type j = 0; j < dofs; j++) {
802 local_ordinal_type blockIndex = i * dofs + j;
803 blockRows[blockIndex].
resize(blockSize);
804 for (
size_t k = 0; k < blockSize; k++) {
807 blockRows[blockIndex][k] = (*Partitioner_)(i, k) * dofs + j;
814 for (local_ordinal_type i = 0; i < NumLocalBlocks_; i++) {
817 "Expected size of all blocks (partitions) to be divisible by MueLu dofs/node.");
818 size_t blockSize = Partitioner_->numRowsInPart(i) / dofs;
822 for (local_ordinal_type j = 0; j < dofs; j++) {
823 local_ordinal_type blockIndex = i * dofs + j;
824 blockRows[blockIndex].
resize(blockSize);
825 for (
size_t k = 0; k < blockSize; k++) {
826 blockRows[blockIndex][k] = (*Partitioner_)(i, k * dofs + j);
833 blockRows.
resize(NumLocalBlocks_);
834 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) {
835 const size_t numRows = Partitioner_->numRowsInPart(i);
836 blockRows[i].
resize(numRows);
838 for (
size_t j = 0; j < numRows; ++j) {
839 blockRows[i][j] = (*Partitioner_)(i, j);
845 Container_->setParameters(List_);
846 Container_->initialize();
849 template <
class MatrixType,
class ContainerType>
850 void BlockRelaxation<MatrixType, ContainerType>::
851 ApplyInverseJacobi(
const MV& X, MV& Y)
const {
852 const size_t NumVectors = X.getNumVectors();
854 MV AY(Y.getMap(), NumVectors);
857 int starting_iteration = 0;
858 if (OverlapLevel_ > 0) {
860 auto WView = W_->getLocalViewHost(Tpetra::Access::ReadOnly);
861 if (ZeroStartingSolution_) {
862 auto XView = X.getLocalViewHost(Tpetra::Access::ReadOnly);
863 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
864 Container_->DoOverlappingJacobi(XView, YView, WView, DampingFactor_, nonsymCombine_);
865 starting_iteration = 1;
867 const scalar_type ONE = STS::one();
868 for (
int j = starting_iteration; j < NumSweeps_; j++) {
870 AY.update(ONE, X, -ONE);
872 auto AYView = AY.getLocalViewHost(Tpetra::Access::ReadOnly);
873 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
874 Container_->DoOverlappingJacobi(AYView, YView, WView, DampingFactor_, nonsymCombine_);
879 if (ZeroStartingSolution_) {
880 auto XView = X.getLocalViewHost(Tpetra::Access::ReadOnly);
881 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
882 Container_->DoJacobi(XView, YView, DampingFactor_);
883 starting_iteration = 1;
885 const scalar_type ONE = STS::one();
886 for (
int j = starting_iteration; j < NumSweeps_; j++) {
888 AY.update(ONE, X, -ONE);
890 auto AYView = AY.getLocalViewHost(Tpetra::Access::ReadOnly);
891 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
892 Container_->DoJacobi(AYView, YView, DampingFactor_);
898 template <
class MatrixType,
class ContainerType>
899 void BlockRelaxation<MatrixType, ContainerType>::
900 ApplyInverseGS(
const MV& X, MV& Y)
const {
903 size_t numVecs = X.getNumVectors();
905 auto XView = X.getLocalViewHost(Tpetra::Access::ReadOnly);
906 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
909 bool deleteY2 =
false;
911 Y2 = ptr(
new MV(Importer_->getTargetMap(), numVecs));
916 for (
int j = 0; j < NumSweeps_; ++j) {
918 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
919 auto Y2View = Y2->getLocalViewHost(Tpetra::Access::ReadWrite);
920 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
923 auto Y2View = Y2->getLocalViewHost(Tpetra::Access::ReadWrite);
924 for (
int j = 0; j < NumSweeps_; ++j) {
925 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
932 template <
class MatrixType,
class ContainerType>
933 void BlockRelaxation<MatrixType, ContainerType>::
934 ApplyInverseSGS(
const MV& X, MV& Y)
const {
938 auto XView = X.getLocalViewHost(Tpetra::Access::ReadOnly);
939 auto YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
942 bool deleteY2 =
false;
944 Y2 = ptr(
new MV(Importer_->getTargetMap(), X.getNumVectors()));
949 for (
int j = 0; j < NumSweeps_; ++j) {
951 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
952 auto Y2View = Y2->getLocalViewHost(Tpetra::Access::ReadWrite);
953 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
956 auto Y2View = Y2->getLocalViewHost(Tpetra::Access::ReadWrite);
957 for (
int j = 0; j < NumSweeps_; ++j) {
958 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
965 template <
class MatrixType,
class ContainerType>
966 void BlockRelaxation<MatrixType, ContainerType>::computeImporter()
const {
973 using Teuchos::rcp_dynamic_cast;
975 if (hasBlockCrsMatrix_) {
976 const RCP<const block_crs_matrix_type> bcrs =
977 rcp_dynamic_cast<
const block_crs_matrix_type>(A_);
978 int bs_ = bcrs->getBlockSize();
979 RCP<const map_type> oldDomainMap = A_->getDomainMap();
980 RCP<const map_type> oldColMap = A_->getColMap();
983 global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
984 global_ordinal_type indexBase = oldColMap->getIndexBase();
985 RCP<const Comm<int>> comm = oldColMap->getComm();
986 ArrayView<const global_ordinal_type> oldColElements = oldColMap->getLocalElementList();
987 Array<global_ordinal_type> newColElements(bs_ * oldColElements.size());
988 for (
int i = 0; i < oldColElements.size(); i++) {
989 for (
int j = 0; j < bs_; j++)
990 newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
992 RCP<map_type> colMap(
new map_type(numGlobalElements, newColElements, indexBase, comm));
994 Importer_ =
rcp(
new import_type(oldDomainMap, colMap));
996 Importer_ = A_->getGraph()->getImporter();
997 if (Importer_.is_null())
998 Importer_ =
rcp(
new import_type(A_->getDomainMap(), A_->getColMap()));
1004 template <
class MatrixType,
class ContainerType>
1008 std::ostringstream out;
1013 out <<
"\"Ifpack2::BlockRelaxation\": {";
1014 if (this->getObjectLabel() !=
"") {
1015 out <<
"Label: \"" << this->getObjectLabel() <<
"\", ";
1020 out <<
"Matrix: null, ";
1024 out <<
"Global matrix dimensions: ["
1025 << A_->getGlobalNumRows() <<
", " << A_->getGlobalNumCols() <<
"], ";
1030 out <<
"\"relaxation: type\": ";
1031 if (PrecType_ == Ifpack2::Details::JACOBI) {
1032 out <<
"Block Jacobi";
1033 }
else if (PrecType_ == Ifpack2::Details::GS) {
1034 out <<
"Block Gauss-Seidel";
1035 }
else if (PrecType_ == Ifpack2::Details::SGS) {
1036 out <<
"Block Symmetric Gauss-Seidel";
1037 }
else if (PrecType_ == Ifpack2::Details::MTSPLITJACOBI) {
1038 out <<
"MT Split Jacobi";
1044 if (hasBlockCrsMatrix_)
1045 out <<
", BlockCrs";
1048 int approx_rows_per_part = A_->getLocalNumRows() / Partitioner_->numLocalParts();
1049 out <<
", blocksize: " << approx_rows_per_part;
1051 out <<
", overlap: " << OverlapLevel_;
1054 <<
"sweeps: " << NumSweeps_ <<
", "
1055 <<
"damping factor: " << DampingFactor_ <<
", ";
1057 std::string containerType = ContainerType::getName();
1058 out <<
"block container: " << ((containerType ==
"Generic") ? containerType_ : containerType);
1059 if (List_.isParameter(
"partitioner: PDE equations"))
1060 out <<
", dofs/node: " << List_.get<
int>(
"partitioner: PDE equations");
1066 template <
class MatrixType,
class ContainerType>
1082 const int myImageID = A_->getComm()->getRank();
1092 if (vl !=
VERB_NONE && myImageID == 0) {
1093 out <<
"Ifpack2::BlockRelaxation<"
1094 << TypeNameTraits<MatrixType>::name() <<
", "
1095 << TypeNameTraits<ContainerType>::name() <<
" >:";
1098 if (this->getObjectLabel() !=
"") {
1099 out <<
"label: \"" << this->getObjectLabel() <<
"\"" << endl;
1101 out <<
"initialized: " << (isInitialized() ?
"true" :
"false") << endl
1102 <<
"computed: " << (isComputed() ?
"true" :
"false") << endl;
1104 std::string precType;
1105 if (PrecType_ == Ifpack2::Details::JACOBI) {
1106 precType =
"Block Jacobi";
1107 }
else if (PrecType_ == Ifpack2::Details::GS) {
1108 precType =
"Block Gauss-Seidel";
1109 }
else if (PrecType_ == Ifpack2::Details::SGS) {
1110 precType =
"Block Symmetric Gauss-Seidel";
1112 out <<
"type: " << precType << endl
1113 <<
"global number of rows: " << A_->getGlobalNumRows() << endl
1114 <<
"global number of columns: " << A_->getGlobalNumCols() << endl
1115 <<
"number of sweeps: " << NumSweeps_ << endl
1116 <<
"damping factor: " << DampingFactor_ << endl
1117 <<
"nonsymmetric overlap combine" << nonsymCombine_ << endl
1118 <<
"backwards mode: "
1119 << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ?
"true" :
"false")
1121 <<
"zero starting solution: "
1122 << (ZeroStartingSolution_ ?
"true" :
"false") << endl;
1139 #ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1141 #define IFPACK2_BLOCKRELAXATION_INSTANT(S, LO, GO, N) \
1142 template class Ifpack2::BlockRelaxation< \
1143 Tpetra::RowMatrix<S, LO, GO, N>, \
1144 Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N>>>;
1146 #endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1148 #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:327
ParameterList & disableRecursiveValidation()
std::string description() const
A one-line description of this object.
Definition: Ifpack2_BlockRelaxation_def.hpp:1007
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:411
T & get(const std::string &name, T def_value)
#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:390
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:26
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:362
int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:384
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void setParameters(const Teuchos::ParameterList ¶ms)
Sets all the parameters for the preconditioner.
Definition: Ifpack2_BlockRelaxation_def.hpp:156
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition: Ifpack2_BlockRelaxation_def.hpp:742
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition: Ifpack2_BlockRelaxation_decl.hpp:49
double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:397
Declaration of a user-defined partitioner in which the user can define a partition of the graph...
Ifpack2::Partitioner:
Definition: Ifpack2_Partitioner.hpp:146
bool isParameter(const std::string &name) const
bool remove(std::string const &name, bool throwIfNotExists=true)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
double getComputeTime() const
Returns the time spent in compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:404
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_BlockRelaxation_def.hpp:416
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_BlockRelaxation_def.hpp:31
virtual ~BlockRelaxation()
Destructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:91
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:349
Partition in which the user can define a nonoverlapping partition of the graph in any way they choose...
Definition: Ifpack2_Details_UserPartitioner_decl.hpp:36
void resize(size_type new_size, const value_type &x=value_type())
int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:378
void initialize()
Initialize.
Definition: Ifpack2_BlockRelaxation_def.hpp:555
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:1068
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:65
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:541
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:62
bool isType(const std::string &name) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:127
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:54
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner's constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:340
void getValidParameters(Teuchos::ParameterList ¶ms)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition: Ifpack2_Parameters.cpp:18
Ifpack2::LinePartitioner: A class to define partitions into a set of lines.
Definition: Ifpack2_LinePartitioner_decl.hpp:45
A class to define linear partitions.
Definition: Ifpack2_LinearPartitioner_decl.hpp:27
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_BlockRelaxation_def.hpp:96
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:433
BlockRelaxation(const Teuchos::RCP< const row_matrix_type > &Matrix)
Constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:58