Ifpack2 Templated Preconditioning Package Version 1.0
Ifpack2_BlockRelaxation_def.hpp
1/*@HEADER
2// ***********************************************************************
3//
4// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5// Copyright (2009) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ***********************************************************************
40//@HEADER
41*/
42
43#ifndef IFPACK2_BLOCKRELAXATION_DEF_HPP
44#define IFPACK2_BLOCKRELAXATION_DEF_HPP
45
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"
53
54namespace Ifpack2 {
55
56template<class MatrixType,class ContainerType>
58setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
59{
60 if (A.getRawPtr () != A_.getRawPtr ()) { // it's a different matrix
61 IsInitialized_ = false;
62 IsComputed_ = false;
63 Partitioner_ = Teuchos::null;
64 W_ = Teuchos::null;
65
66 if (! A.is_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();
71 }
72 if (! Container_.is_null ()) {
73 //This just frees up the container's memory.
74 //The container will be fully re-constructed during initialize().
75 Container_->clearBlocks();
76 }
77 NumLocalBlocks_ = 0;
78
79 A_ = A;
80 computeImporter();
81 }
82}
83
84template<class MatrixType,class ContainerType>
86BlockRelaxation (const Teuchos::RCP<const row_matrix_type>& A)
87:
88 Container_ (Teuchos::null),
89 Partitioner_ (Teuchos::null),
90 PartitionerType_ ("linear"),
91 NumSweeps_ (1),
92 NumLocalBlocks_(0),
93 containerType_ ("TriDi"),
94 PrecType_ (Ifpack2::Details::JACOBI),
95 ZeroStartingSolution_ (true),
96 hasBlockCrsMatrix_ (false),
97 DoBackwardGS_ (false),
98 OverlapLevel_ (0),
99 DampingFactor_ (STS::one ()),
100 IsInitialized_ (false),
101 IsComputed_ (false),
102 NumInitialize_ (0),
103 NumCompute_ (0),
104 NumApply_ (0),
105 InitializeTime_ (0.0),
106 ComputeTime_ (0.0),
107 NumLocalRows_ (0),
108 NumGlobalRows_ (0),
109 NumGlobalNonzeros_ (0),
110 W_ (Teuchos::null),
111 Importer_ (Teuchos::null)
112{
113 setMatrix(A);
114}
115
116template<class MatrixType,class ContainerType>
119{}
120
121template<class MatrixType,class ContainerType>
122Teuchos::RCP<const Teuchos::ParameterList>
124getValidParameters () const
125{
126 Teuchos::RCP<Teuchos::ParameterList> validParams = Teuchos::parameterList ("Ifpack2::BlockRelaxation");
127
128 validParams->set("relaxation: container", "TriDi");
129 validParams->set("relaxation: backward mode",false);
130 validParams->set("relaxation: type", "Jacobi");
131 validParams->set("relaxation: sweeps", 1);
132 validParams->set("relaxation: damping factor", STS::one());
133 validParams->set("relaxation: zero starting solution", true);
134 validParams->set("block relaxation: decouple dofs", false);
135 validParams->set("schwarz: compute condest", false); // mfh 24 Mar 2015: for backwards compatibility ONLY
136 validParams->set("schwarz: combine mode", "ZERO"); // use string mode for this
137 validParams->set("schwarz: use reordering", true);
138 validParams->set("schwarz: filter singletons", false);
139 validParams->set("schwarz: overlap level", 0);
140 validParams->set("partitioner: type", "greedy");
141 validParams->set("partitioner: local parts", 1);
142 validParams->set("partitioner: overlap", 0);
143 Teuchos::Array<Teuchos::ArrayRCP<int>> tmp0;
144 validParams->set("partitioner: parts", tmp0);
145 validParams->set("partitioner: maintain sparsity", false);
146 validParams->set("fact: ilut level-of-fill", 1.0);
147 validParams->set("fact: absolute threshold", 0.0);
148 validParams->set("fact: relative threshold", 1.0);
149 validParams->set("fact: relax value", 0.0);
150
151 Teuchos::ParameterList dummyList;
152 validParams->set("Amesos2",dummyList);
153 validParams->sublist("Amesos2").disableRecursiveValidation();
154 validParams->set("Amesos2 solver name", "KLU2");
155
156 Teuchos::ArrayRCP<int> tmp;
157 validParams->set("partitioner: map", tmp);
158
159 validParams->set("partitioner: line detection threshold", 0.0);
160 validParams->set("partitioner: PDE equations", 1);
161 Teuchos::RCP<Tpetra::MultiVector<double,
162 typename MatrixType::local_ordinal_type,
163 typename MatrixType::global_ordinal_type,
164 typename MatrixType::node_type> > dummy;
165 validParams->set("partitioner: coordinates",dummy);
166
167 return validParams;
168}
169
170template<class MatrixType,class ContainerType>
171void
173setParameters (const Teuchos::ParameterList& List)
174{
175 // Note that the validation process does not change List.
176 Teuchos::RCP<const Teuchos::ParameterList> validparams;
177 validparams = this->getValidParameters();
178 List.validateParameters (*validparams);
179
180 if (List.isParameter ("relaxation: container")) {
181 // If the container type isn't a string, this will throw, but it
182 // rightfully should.
183
184 // If its value does not match the currently registered Container types,
185 // the ContainerFactory will throw with an informative message.
186 containerType_ = List.get<std::string> ("relaxation: container");
187 // Intercept more human-readable aliases for the *TriDi containers
188 if(containerType_ == "Tridiagonal") {
189 containerType_ = "TriDi";
190 }
191 if(containerType_ == "Block Tridiagonal") {
192 containerType_ = "BlockTriDi";
193 }
194 }
195
196 if (List.isParameter ("relaxation: type")) {
197 const std::string relaxType = List.get<std::string> ("relaxation: type");
198
199 if (relaxType == "Jacobi") {
200 PrecType_ = Ifpack2::Details::JACOBI;
201 }
202 else if (relaxType == "MT Split Jacobi") {
203 PrecType_ = Ifpack2::Details::MTSPLITJACOBI;
204 }
205 else if (relaxType == "Gauss-Seidel") {
206 PrecType_ = Ifpack2::Details::GS;
207 }
208 else if (relaxType == "Symmetric Gauss-Seidel") {
209 PrecType_ = Ifpack2::Details::SGS;
210 }
211 else {
212 TEUCHOS_TEST_FOR_EXCEPTION
213 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
214 "setParameters: Invalid parameter value \"" << relaxType
215 << "\" for parameter \"relaxation: type\".");
216 }
217 }
218
219 if (List.isParameter ("relaxation: sweeps")) {
220 NumSweeps_ = List.get<int> ("relaxation: sweeps");
221 }
222
223 // Users may specify this as a floating-point literal. In that
224 // case, it may have type double, even if scalar_type is something
225 // else. We take care to try the various types that make sense.
226 if (List.isParameter ("relaxation: damping factor")) {
227 if (List.isType<double> ("relaxation: damping factor")) {
228 const double dampingFactor =
229 List.get<double> ("relaxation: damping factor");
230 DampingFactor_ = static_cast<scalar_type> (dampingFactor);
231 }
232 else if (List.isType<scalar_type> ("relaxation: damping factor")) {
233 DampingFactor_ = List.get<scalar_type> ("relaxation: damping factor");
234 }
235 else if (List.isType<magnitude_type> ("relaxation: damping factor")) {
236 const magnitude_type dampingFactor =
237 List.get<magnitude_type> ("relaxation: damping factor");
238 DampingFactor_ = static_cast<scalar_type> (dampingFactor);
239 }
240 else {
241 TEUCHOS_TEST_FOR_EXCEPTION
242 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
243 "setParameters: Parameter \"relaxation: damping factor\" "
244 "has the wrong type.");
245 }
246 }
247
248 if (List.isParameter ("relaxation: zero starting solution")) {
249 ZeroStartingSolution_ = List.get<bool> ("relaxation: zero starting solution");
250 }
251
252 if (List.isParameter ("relaxation: backward mode")) {
253 DoBackwardGS_ = List.get<bool> ("relaxation: backward mode");
254 }
255
256 if (List.isParameter ("partitioner: type")) {
257 PartitionerType_ = List.get<std::string> ("partitioner: type");
258 }
259
260 // Users may specify this as an int literal, so we need to allow
261 // both int and local_ordinal_type (not necessarily same as int).
262 if (List.isParameter ("partitioner: local parts")) {
263 if (List.isType<local_ordinal_type> ("partitioner: local parts")) {
264 NumLocalBlocks_ = List.get<local_ordinal_type> ("partitioner: local parts");
265 }
266 else if (! std::is_same<int, local_ordinal_type>::value &&
267 List.isType<int> ("partitioner: local parts")) {
268 NumLocalBlocks_ = List.get<int> ("partitioner: local parts");
269 }
270 else {
271 TEUCHOS_TEST_FOR_EXCEPTION
272 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
273 "setParameters: Parameter \"partitioner: local parts\" "
274 "has the wrong type.");
275 }
276 }
277
278 if (List.isParameter ("partitioner: overlap level")) {
279 if (List.isType<int> ("partitioner: overlap level")) {
280 OverlapLevel_ = List.get<int> ("partitioner: overlap level");
281 }
282 else if (! std::is_same<int, local_ordinal_type>::value &&
283 List.isType<local_ordinal_type> ("partitioner: overlap level")) {
284 OverlapLevel_ = List.get<local_ordinal_type> ("partitioner: overlap level");
285 }
286 else {
287 TEUCHOS_TEST_FOR_EXCEPTION
288 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
289 "setParameters: Parameter \"partitioner: overlap level\" "
290 "has the wrong type.");
291 }
292 }
293
294 std::string defaultContainer = "TriDi";
295 if(std::is_same<ContainerType, Container<MatrixType> >::value)
296 {
297 //Generic container template parameter, container type specified in List
298 Ifpack2::getParameter(List, "relaxation: container", defaultContainer);
299 }
300 // check parameters
301 if (PrecType_ != Ifpack2::Details::JACOBI) {
302 OverlapLevel_ = 0;
303 }
304 if (NumLocalBlocks_ < static_cast<local_ordinal_type> (0)) {
305 NumLocalBlocks_ = A_->getNodeNumRows() / (-NumLocalBlocks_);
306 }
307
308 decouple_ = false;
309 if(List.isParameter("block relaxation: decouple dofs"))
310 decouple_ = List.get<bool>("block relaxation: decouple dofs");
311
312 // other checks are performed in Partitioner_
313
314 // NTS: Sanity check to be removed at a later date when Backward mode is enabled
315 TEUCHOS_TEST_FOR_EXCEPTION(
316 DoBackwardGS_, std::runtime_error,
317 "Ifpack2::BlockRelaxation:setParameters: Setting the \"relaxation: "
318 "backward mode\" parameter to true is not yet supported.");
319
320 // copy the list as each subblock's constructor will
321 // require it later
322 List_ = List;
323}
324
325template<class MatrixType,class ContainerType>
326Teuchos::RCP<const Teuchos::Comm<int> >
328{
329 TEUCHOS_TEST_FOR_EXCEPTION
330 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::getComm: "
331 "The matrix is null. You must call setMatrix() with a nonnull matrix "
332 "before you may call this method.");
333 return A_->getComm ();
334}
335
336template<class MatrixType,class ContainerType>
337Teuchos::RCP<const Tpetra::RowMatrix<typename MatrixType::scalar_type,
338 typename MatrixType::local_ordinal_type,
339 typename MatrixType::global_ordinal_type,
340 typename MatrixType::node_type> >
342 return A_;
343}
344
345template<class MatrixType,class ContainerType>
346Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
347 typename MatrixType::global_ordinal_type,
348 typename MatrixType::node_type> >
350getDomainMap () const
351{
352 TEUCHOS_TEST_FOR_EXCEPTION
353 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
354 "getDomainMap: The matrix is null. You must call setMatrix() with a "
355 "nonnull matrix before you may call this method.");
356 return A_->getDomainMap ();
357}
358
359template<class MatrixType,class ContainerType>
360Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
361 typename MatrixType::global_ordinal_type,
362 typename MatrixType::node_type> >
364getRangeMap () const
365{
366 TEUCHOS_TEST_FOR_EXCEPTION
367 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
368 "getRangeMap: The matrix is null. You must call setMatrix() with a "
369 "nonnull matrix before you may call this method.");
370 return A_->getRangeMap ();
371}
372
373template<class MatrixType,class ContainerType>
374bool
376hasTransposeApply () const {
377 return true;
378}
379
380template<class MatrixType,class ContainerType>
381int
383getNumInitialize () const {
384 return NumInitialize_;
385}
386
387template<class MatrixType,class ContainerType>
388int
390getNumCompute () const
391{
392 return NumCompute_;
393}
394
395template<class MatrixType,class ContainerType>
396int
398getNumApply () const
399{
400 return NumApply_;
401}
402
403template<class MatrixType,class ContainerType>
404double
406getInitializeTime () const
407{
408 return InitializeTime_;
409}
410
411template<class MatrixType,class ContainerType>
412double
414getComputeTime () const
415{
416 return ComputeTime_;
417}
418
419template<class MatrixType,class ContainerType>
420double
422getApplyTime () const
423{
424 return ApplyTime_;
425}
426
427
428template<class MatrixType,class ContainerType>
430 TEUCHOS_TEST_FOR_EXCEPTION(
431 A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::getNodeSmootherComplexity: "
432 "The input matrix A is null. Please call setMatrix() with a nonnull "
433 "input matrix, then call compute(), before calling this method.");
434 // Relaxation methods cost roughly one apply + one block-diagonal inverse per iteration
435 // NOTE: This approximates all blocks as dense, which may overstate the cost if you have a sparse (or banded) container.
436 size_t block_nnz = 0;
437 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i)
438 block_nnz += Partitioner_->numRowsInPart(i) *Partitioner_->numRowsInPart(i);
439
440 return block_nnz + A_->getNodeNumEntries();
441}
442
443template<class MatrixType,class ContainerType>
444void
446apply (const Tpetra::MultiVector<typename MatrixType::scalar_type,
447 typename MatrixType::local_ordinal_type,
448 typename MatrixType::global_ordinal_type,
449 typename MatrixType::node_type>& X,
450 Tpetra::MultiVector<typename MatrixType::scalar_type,
451 typename MatrixType::local_ordinal_type,
452 typename MatrixType::global_ordinal_type,
453 typename MatrixType::node_type>& Y,
454 Teuchos::ETransp mode,
455 scalar_type alpha,
456 scalar_type beta) const
457{
458 TEUCHOS_TEST_FOR_EXCEPTION
459 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::apply: "
460 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
461 "then call initialize() and compute() (in that order), before you may "
462 "call this method.");
463 TEUCHOS_TEST_FOR_EXCEPTION(
464 ! isComputed (), std::runtime_error, "Ifpack2::BlockRelaxation::apply: "
465 "isComputed() must be true prior to calling apply.");
466 TEUCHOS_TEST_FOR_EXCEPTION(
467 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
468 "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = "
469 << X.getNumVectors () << " != Y.getNumVectors() = "
470 << Y.getNumVectors () << ".");
471 TEUCHOS_TEST_FOR_EXCEPTION(
472 mode != Teuchos::NO_TRANS, std::logic_error, "Ifpack2::BlockRelaxation::"
473 "apply: This method currently only implements the case mode == Teuchos::"
474 "NO_TRANS. Transposed modes are not currently supported.");
475 TEUCHOS_TEST_FOR_EXCEPTION(
476 alpha != Teuchos::ScalarTraits<scalar_type>::one (), std::logic_error,
477 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
478 "the case alpha == 1. You specified alpha = " << alpha << ".");
479 TEUCHOS_TEST_FOR_EXCEPTION(
480 beta != Teuchos::ScalarTraits<scalar_type>::zero (), std::logic_error,
481 "Ifpack2::BlockRelaxation::apply: This method currently only implements "
482 "the case beta == 0. You specified beta = " << beta << ".");
483
484 const std::string timerName ("Ifpack2::BlockRelaxation::apply");
485 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
486 if (timer.is_null ()) {
487 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
488 }
489
490 double startTime = timer->wallTime();
491
492 {
493 Teuchos::TimeMonitor timeMon (*timer);
494
495 // If X and Y are pointing to the same memory location,
496 // we need to create an auxiliary vector, Xcopy
497 Teuchos::RCP<const MV> X_copy;
498 {
499 if (X.aliases(Y)) {
500 X_copy = rcp (new MV (X, Teuchos::Copy));
501 } else {
502 X_copy = rcpFromRef (X);
503 }
504 }
505
506 if (ZeroStartingSolution_) {
507 Y.putScalar (STS::zero ());
508 }
509
510 // Flops are updated in each of the following.
511 switch (PrecType_) {
512 case Ifpack2::Details::JACOBI:
513 ApplyInverseJacobi(*X_copy,Y);
514 break;
515 case Ifpack2::Details::GS:
516 ApplyInverseGS(*X_copy,Y);
517 break;
518 case Ifpack2::Details::SGS:
519 ApplyInverseSGS(*X_copy,Y);
520 break;
521 case Ifpack2::Details::MTSPLITJACOBI:
522 //note: for this method, the container is always BlockTriDi
523 Container_->applyInverseJacobi(*X_copy, Y, DampingFactor_, ZeroStartingSolution_, NumSweeps_);
524 break;
525 default:
526 TEUCHOS_TEST_FOR_EXCEPTION
527 (true, std::logic_error, "Ifpack2::BlockRelaxation::apply: Invalid "
528 "PrecType_ enum value " << PrecType_ << ". Valid values are Ifpack2::"
529 "Details::JACOBI = " << Ifpack2::Details::JACOBI << ", Ifpack2::Details"
530 "::GS = " << Ifpack2::Details::GS << ", and Ifpack2::Details::SGS = "
531 << Ifpack2::Details::SGS << ". Please report this bug to the Ifpack2 "
532 "developers.");
533 }
534 }
535
536 ApplyTime_ += (timer->wallTime() - startTime);
537 ++NumApply_;
538}
539
540template<class MatrixType,class ContainerType>
541void
543applyMat (const Tpetra::MultiVector<typename MatrixType::scalar_type,
544 typename MatrixType::local_ordinal_type,
545 typename MatrixType::global_ordinal_type,
546 typename MatrixType::node_type>& X,
547 Tpetra::MultiVector<typename MatrixType::scalar_type,
548 typename MatrixType::local_ordinal_type,
549 typename MatrixType::global_ordinal_type,
550 typename MatrixType::node_type>& Y,
551 Teuchos::ETransp mode) const
552{
553 A_->apply (X, Y, mode);
554}
555
556template<class MatrixType,class ContainerType>
557void
559initialize ()
560{
561 using Teuchos::rcp;
562 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
563 row_graph_type;
564
565 TEUCHOS_TEST_FOR_EXCEPTION
566 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::initialize: "
567 "The matrix is null. You must call setMatrix() with a nonnull matrix "
568 "before you may call this method.");
569
570 Teuchos::RCP<Teuchos::Time> timer =
571 Teuchos::TimeMonitor::getNewCounter ("Ifpack2::BlockRelaxation::initialize");
572 double startTime = timer->wallTime();
573
574 { // Timing of initialize starts here
575 Teuchos::TimeMonitor timeMon (*timer);
576 IsInitialized_ = false;
577
578 // Check whether we have a BlockCrsMatrix
579 Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
580 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A_);
581 hasBlockCrsMatrix_ = !A_bcrs.is_null();
582 if (A_bcrs.is_null ()) {
583 hasBlockCrsMatrix_ = false;
584 }
585 else {
586 hasBlockCrsMatrix_ = true;
587 }
588
589 NumLocalRows_ = A_->getNodeNumRows ();
590 NumGlobalRows_ = A_->getGlobalNumRows ();
591 NumGlobalNonzeros_ = A_->getGlobalNumEntries ();
592
593 // NTS: Will need to add support for Zoltan2 partitions later Also,
594 // will need a better way of handling the Graph typing issue.
595 // Especially with ordinal types w.r.t the container.
596 Partitioner_ = Teuchos::null;
597
598 if (PartitionerType_ == "linear") {
599 Partitioner_ =
600 rcp (new Ifpack2::LinearPartitioner<row_graph_type> (A_->getGraph ()));
601 } else if (PartitionerType_ == "line") {
602 Partitioner_ =
604 } else if (PartitionerType_ == "user") {
605 Partitioner_ =
606 rcp (new Ifpack2::Details::UserPartitioner<row_graph_type> (A_->getGraph () ) );
607 } else {
608 // We should have checked for this in setParameters(), so it's a
609 // logic_error, not an invalid_argument or runtime_error.
610 TEUCHOS_TEST_FOR_EXCEPTION
611 (true, std::logic_error, "Ifpack2::BlockRelaxation::initialize: Unknown "
612 "partitioner type " << PartitionerType_ << ". Valid values are "
613 "\"linear\", \"line\", and \"user\".");
614 }
615
616 // need to partition the graph of A
617 Partitioner_->setParameters (List_);
618 Partitioner_->compute ();
619
620 // get actual number of partitions
621 NumLocalBlocks_ = Partitioner_->numLocalParts ();
622
623 // Note: Unlike Ifpack, we'll punt on computing W_ until compute(), which is where
624 // we assume that the type of relaxation has been chosen.
625
626 if (A_->getComm()->getSize() != 1) {
627 IsParallel_ = true;
628 } else {
629 IsParallel_ = false;
630 }
631
632 // We should have checked for this in setParameters(), so it's a
633 // logic_error, not an invalid_argument or runtime_error.
634 TEUCHOS_TEST_FOR_EXCEPTION
635 (NumSweeps_ < 0, std::logic_error, "Ifpack2::BlockRelaxation::initialize: "
636 "NumSweeps_ = " << NumSweeps_ << " < 0.");
637
638 // Extract the submatrices
639 ExtractSubmatricesStructure ();
640
641 // Compute the weight vector if we're doing overlapped Jacobi (and
642 // only if we're doing overlapped Jacobi).
643 if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) {
644 TEUCHOS_TEST_FOR_EXCEPTION
645 (hasBlockCrsMatrix_, std::runtime_error,
646 "Ifpack2::BlockRelaxation::initialize: "
647 "We do not support overlapped Jacobi yet for Tpetra::BlockCrsMatrix. Sorry!");
648
649 // weight of each vertex
650 W_ = rcp (new vector_type (A_->getRowMap ()));
651 W_->putScalar (STS::zero ());
652 Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0);
653
654 for (local_ordinal_type i = 0 ; i < NumLocalBlocks_ ; ++i) {
655 for (size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) {
656 local_ordinal_type LID = (*Partitioner_)(i,j);
657 w_ptr[LID] += STS::one();
658 }
659 }
660 W_->reciprocal (*W_);
661 }
662 } // timing of initialize stops here
663
664 InitializeTime_ += (timer->wallTime() - startTime);
665 ++NumInitialize_;
666 IsInitialized_ = true;
667}
668
669
670template<class MatrixType,class ContainerType>
671void
673compute ()
674{
675 using Teuchos::rcp;
676
677 TEUCHOS_TEST_FOR_EXCEPTION
678 (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::compute: "
679 "The matrix is null. You must call setMatrix() with a nonnull matrix, "
680 "then call initialize(), before you may call this method.");
681
682 if (! isInitialized ()) {
683 initialize ();
684 }
685
686 Teuchos::RCP<Teuchos::Time> timer =
687 Teuchos::TimeMonitor::getNewCounter ("Ifpack2::BlockRelaxation::compute");
688
689 double startTime = timer->wallTime();
690
691 {
692 Teuchos::TimeMonitor timeMon (*timer);
693
694 // reset values
695 IsComputed_ = false;
696
697 Container_->compute(); // compute each block matrix
698 }
699
700 ComputeTime_ += (timer->wallTime() - startTime);
701 ++NumCompute_;
702 IsComputed_ = true;
703}
704
705template<class MatrixType, class ContainerType>
706void
709{
710 TEUCHOS_TEST_FOR_EXCEPTION
711 (Partitioner_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
712 "ExtractSubmatricesStructure: Partitioner object is null.");
713
714 std::string containerType = ContainerType::getName ();
715 if (containerType == "Generic") {
716 // ContainerType is Container<row_matrix_type>. Thus, we need to
717 // get the container name from the parameter list. We use "TriDi"
718 // as the default container type.
719 containerType = containerType_;
720 }
721 //Whether the Container will define blocks (partitions)
722 //in terms of individual DOFs, and not by nodes (blocks).
723 bool pointIndexed = decouple_ && hasBlockCrsMatrix_;
724 Teuchos::Array<Teuchos::Array<local_ordinal_type> > blockRows;
725 if(decouple_)
726 {
727 //dofs [per node] is how many blocks each partition will be split into
728 auto A_bcrs = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(A_);
729 local_ordinal_type dofs = hasBlockCrsMatrix_ ?
730 A_bcrs->getBlockSize() : List_.get<int>("partitioner: PDE equations");
731 blockRows.resize(NumLocalBlocks_ * dofs);
732 if(hasBlockCrsMatrix_)
733 {
734 for(local_ordinal_type i = 0; i < NumLocalBlocks_; i++)
735 {
736 size_t blockSize = Partitioner_->numRowsInPart(i);
737 //block i will be split into j different blocks,
738 //each corresponding to a different dof
739 for(local_ordinal_type j = 0; j < dofs; j++)
740 {
741 local_ordinal_type blockIndex = i * dofs + j;
742 blockRows[blockIndex].resize(blockSize);
743 for(size_t k = 0; k < blockSize; k++)
744 {
745 //here, the row and dof are combined to get the point index
746 //(what the row would be if A were a CrsMatrix)
747 blockRows[blockIndex][k] = (*Partitioner_)(i, k) * dofs + j;
748 }
749 }
750 }
751 }
752 else
753 {
754 //Rows in each partition are distributed round-robin to the blocks -
755 //that's how MueLu stores DOFs in a non-block matrix
756 for(local_ordinal_type i = 0; i < NumLocalBlocks_; i++)
757 {
758 //#ifdef HAVE_IFPACK2_DEBUG
759 TEUCHOS_TEST_FOR_EXCEPTION(Partitioner_->numRowsInPart(i) % dofs != 0, std::logic_error,
760 "Expected size of all blocks (partitions) to be divisible by MueLu dofs/node.");
761 size_t blockSize = Partitioner_->numRowsInPart(i) / dofs;
762 //#endif
763 //block i will be split into j different blocks,
764 //each corresponding to a different dof
765 for(local_ordinal_type j = 0; j < dofs; j++)
766 {
767 local_ordinal_type blockIndex = i * dofs + j;
768 blockRows[blockIndex].resize(blockSize);
769 for(size_t k = 0; k < blockSize; k++)
770 {
771 blockRows[blockIndex][k] = (*Partitioner_)(i, k * dofs + j);
772 }
773 }
774 }
775 }
776 }
777 else
778 {
779 //No decoupling - just get the rows directly from Partitioner
780 blockRows.resize(NumLocalBlocks_);
781 for(local_ordinal_type i = 0; i < NumLocalBlocks_; ++i)
782 {
783 const size_t numRows = Partitioner_->numRowsInPart (i);
784 blockRows[i].resize(numRows);
785 // Extract a list of the indices of each partitioner row.
786 for (size_t j = 0; j < numRows; ++j)
787 {
788 blockRows[i][j] = (*Partitioner_) (i,j);
789 }
790 }
791 }
792 //right before constructing the
793 Container_ = ContainerFactory<MatrixType>::build(containerType, A_, blockRows, Importer_, pointIndexed);
794 Container_->setParameters(List_);
795 Container_->initialize();
796}
797
798template<class MatrixType,class ContainerType>
799void
800BlockRelaxation<MatrixType,ContainerType>::
801ApplyInverseJacobi (const MV& X, MV& Y) const
802{
803 const size_t NumVectors = X.getNumVectors ();
804
805 MV AY (Y.getMap (), NumVectors);
806
807 // Initial matvec not needed
808 int starting_iteration = 0;
809 if (OverlapLevel_ > 0)
810 {
811 //Overlapping jacobi, with view of W_
812 auto WView = W_->getLocalViewHost (Tpetra::Access::ReadOnly);
813 if(ZeroStartingSolution_) {
814 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
815 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
816 Container_->DoOverlappingJacobi(XView, YView, WView, DampingFactor_);
817 starting_iteration = 1;
818 }
819 const scalar_type ONE = STS::one();
820 for(int j = starting_iteration; j < NumSweeps_; j++)
821 {
822 applyMat (Y, AY);
823 AY.update (ONE, X, -ONE);
824 {
825 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
826 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
827 Container_->DoOverlappingJacobi (AYView, YView, WView, DampingFactor_);
828 }
829 }
830 }
831 else
832 {
833 //Non-overlapping
834 if(ZeroStartingSolution_)
835 {
836 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
837 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
838 Container_->DoJacobi (XView, YView, DampingFactor_);
839 starting_iteration = 1;
840 }
841 const scalar_type ONE = STS::one();
842 for(int j = starting_iteration; j < NumSweeps_; j++)
843 {
844 applyMat (Y, AY);
845 AY.update (ONE, X, -ONE);
846 {
847 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
848 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
849 Container_->DoJacobi (AYView, YView, DampingFactor_);
850 }
851 }
852 }
853}
854
855template<class MatrixType,class ContainerType>
856void
857BlockRelaxation<MatrixType,ContainerType>::
858ApplyInverseGS (const MV& X, MV& Y) const
859{
860 using Teuchos::Ptr;
861 using Teuchos::ptr;
862 size_t numVecs = X.getNumVectors();
863 //Get view of X (is never modified in this function)
864 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
865 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
866 //Pre-import Y, if parallel
867 Ptr<MV> Y2;
868 bool deleteY2 = false;
869 if(IsParallel_)
870 {
871 Y2 = ptr(new MV(Importer_->getTargetMap(), numVecs));
872 deleteY2 = true;
873 }
874 else
875 Y2 = ptr(&Y);
876 if(IsParallel_)
877 {
878 for(int j = 0; j < NumSweeps_; ++j)
879 {
880 //do import once per sweep
881 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
882 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
883 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
884 }
885 }
886 else
887 {
888 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
889 for(int j = 0; j < NumSweeps_; ++j)
890 {
891 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
892 }
893 }
894 if(deleteY2)
895 delete Y2.get();
896}
897
898template<class MatrixType,class ContainerType>
899void
900BlockRelaxation<MatrixType,ContainerType>::
901ApplyInverseSGS (const MV& X, MV& Y) const
902{
903 using Teuchos::Ptr;
904 using Teuchos::ptr;
905 //Get view of X (is never modified in this function)
906 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
907 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
908 //Pre-import Y, if parallel
909 Ptr<MV> Y2;
910 bool deleteY2 = false;
911 if(IsParallel_)
912 {
913 Y2 = ptr(new MV(Importer_->getTargetMap(), X.getNumVectors()));
914 deleteY2 = true;
915 }
916 else
917 Y2 = ptr(&Y);
918 if(IsParallel_)
919 {
920 for(int j = 0; j < NumSweeps_; ++j)
921 {
922 //do import once per sweep
923 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
924 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
925 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
926 }
927 }
928 else
929 {
930 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
931 for(int j = 0; j < NumSweeps_; ++j)
932 {
933 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
934 }
935 }
936 if(deleteY2)
937 delete Y2.get();
938}
939
940template<class MatrixType,class ContainerType>
941void BlockRelaxation<MatrixType,ContainerType>::computeImporter () const
942{
943 using Teuchos::RCP;
944 using Teuchos::rcp;
945 using Teuchos::Ptr;
946 using Teuchos::ArrayView;
947 using Teuchos::Array;
948 using Teuchos::Comm;
949 using Teuchos::rcp_dynamic_cast;
950 if(IsParallel_)
951 {
952 if(hasBlockCrsMatrix_)
953 {
954 const RCP<const block_crs_matrix_type> bcrs =
955 rcp_dynamic_cast<const block_crs_matrix_type>(A_);
956 int bs_ = bcrs->getBlockSize();
957 RCP<const map_type> oldDomainMap = A_->getDomainMap();
958 RCP<const map_type> oldColMap = A_->getColMap();
959 // Because A is a block CRS matrix, import will not do what you think it does
960 // We have to construct the correct maps for it
961 global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
962 global_ordinal_type indexBase = oldColMap->getIndexBase();
963 RCP<const Comm<int>> comm = oldColMap->getComm();
964 ArrayView<const global_ordinal_type> oldColElements = oldColMap->getNodeElementList();
965 Array<global_ordinal_type> newColElements(bs_ * oldColElements.size());
966 for(int i = 0; i < oldColElements.size(); i++)
967 {
968 for(int j = 0; j < bs_; j++)
969 newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
970 }
971 RCP<map_type> colMap(new map_type(numGlobalElements, newColElements, indexBase, comm));
972 // Create the importer
973 Importer_ = rcp(new import_type(oldDomainMap, colMap));
974 }
975 else if(!A_.is_null())
976 {
977 Importer_ = A_->getGraph()->getImporter();
978 if(Importer_.is_null())
979 Importer_ = rcp(new import_type(A_->getDomainMap(), A_->getColMap()));
980 }
981 }
982 //otherwise, Importer_ is not needed and is left NULL
983}
984
985template<class MatrixType, class ContainerType>
986std::string
988description () const
989{
990 std::ostringstream out;
991
992 // Output is a valid YAML dictionary in flow style. If you don't
993 // like everything on a single line, you should call describe()
994 // instead.
995 out << "\"Ifpack2::BlockRelaxation\": {";
996 if (this->getObjectLabel () != "") {
997 out << "Label: \"" << this->getObjectLabel () << "\", ";
998 }
999 // out << "Initialized: " << (isInitialized () ? "true" : "false") << ", ";
1000 // out << "Computed: " << (isComputed () ? "true" : "false") << ", ";
1001 if (A_.is_null ()) {
1002 out << "Matrix: null, ";
1003 }
1004 else {
1005 // out << "Matrix: not null"
1006 // << ", Global matrix dimensions: ["
1007 out << "Global matrix dimensions: ["
1008 << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "], ";
1009 }
1010
1011 // It's useful to print this instance's relaxation method. If you
1012 // want more info than that, call describe() instead.
1013 out << "\"relaxation: type\": ";
1014 if (PrecType_ == Ifpack2::Details::JACOBI) {
1015 out << "Block Jacobi";
1016 } else if (PrecType_ == Ifpack2::Details::GS) {
1017 out << "Block Gauss-Seidel";
1018 } else if (PrecType_ == Ifpack2::Details::SGS) {
1019 out << "Block Symmetric Gauss-Seidel";
1020 } else if (PrecType_ == Ifpack2::Details::MTSPLITJACOBI) {
1021 out << "MT Split Jacobi";
1022 } else {
1023 out << "INVALID";
1024 }
1025 // Print the approximate # rows per part
1026 int approx_rows_per_part = A_->getNodeNumRows()/Partitioner_->numLocalParts();
1027 out <<", blocksize: "<<approx_rows_per_part;
1028
1029 out << ", overlap: " << OverlapLevel_;
1030
1031 out << ", " << "sweeps: " << NumSweeps_ << ", "
1032 << "damping factor: " << DampingFactor_ << ", ";
1033
1034 std::string containerType = ContainerType::getName();
1035 out << "block container: " << ((containerType == "Generic") ? containerType_ : containerType);
1036 if(List_.isParameter("partitioner: PDE equations"))
1037 out << ", dofs/node: "<<List_.get<int>("partitioner: PDE equations");
1038
1039
1040 out << "}";
1041 return out.str();
1042}
1043
1044template<class MatrixType,class ContainerType>
1045void
1047describe (Teuchos::FancyOStream& out,
1048 const Teuchos::EVerbosityLevel verbLevel) const
1049{
1050 using std::endl;
1051 using std::setw;
1052 using Teuchos::TypeNameTraits;
1053 using Teuchos::VERB_DEFAULT;
1054 using Teuchos::VERB_NONE;
1055 using Teuchos::VERB_LOW;
1056 using Teuchos::VERB_MEDIUM;
1057 using Teuchos::VERB_HIGH;
1058 using Teuchos::VERB_EXTREME;
1059
1060 Teuchos::EVerbosityLevel vl = verbLevel;
1061 if (vl == VERB_DEFAULT) vl = VERB_LOW;
1062 const int myImageID = A_->getComm()->getRank();
1063
1064 // Convention requires that describe() begin with a tab.
1065 Teuchos::OSTab tab (out);
1066
1067 // none: print nothing
1068 // low: print O(1) info from node 0
1069 // medium:
1070 // high:
1071 // extreme:
1072 if (vl != VERB_NONE && myImageID == 0) {
1073 out << "Ifpack2::BlockRelaxation<"
1074 << TypeNameTraits<MatrixType>::name () << ", "
1075 << TypeNameTraits<ContainerType>::name () << " >:";
1076 Teuchos::OSTab tab1 (out);
1077
1078 if (this->getObjectLabel () != "") {
1079 out << "label: \"" << this->getObjectLabel () << "\"" << endl;
1080 }
1081 out << "initialized: " << (isInitialized () ? "true" : "false") << endl
1082 << "computed: " << (isComputed () ? "true" : "false") << endl;
1083
1084 std::string precType;
1085 if (PrecType_ == Ifpack2::Details::JACOBI) {
1086 precType = "Block Jacobi";
1087 } else if (PrecType_ == Ifpack2::Details::GS) {
1088 precType = "Block Gauss-Seidel";
1089 } else if (PrecType_ == Ifpack2::Details::SGS) {
1090 precType = "Block Symmetric Gauss-Seidel";
1091 }
1092 out << "type: " << precType << endl
1093 << "global number of rows: " << A_->getGlobalNumRows () << endl
1094 << "global number of columns: " << A_->getGlobalNumCols () << endl
1095 << "number of sweeps: " << NumSweeps_ << endl
1096 << "damping factor: " << DampingFactor_ << endl
1097 << "backwards mode: "
1098 << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ? "true" : "false")
1099 << endl
1100 << "zero starting solution: "
1101 << (ZeroStartingSolution_ ? "true" : "false") << endl;
1102 }
1103}
1104
1105}//namespace Ifpack2
1106
1107
1108// Macro that does explicit template instantiation (ETI) for
1109// Ifpack2::BlockRelaxation. S, LO, GO, N correspond to the four
1110// template parameters of Ifpack2::Preconditioner and
1111// Tpetra::RowMatrix.
1112//
1113// We only instantiate for MatrixType = Tpetra::RowMatrix. There's no
1114// need to instantiate for Tpetra::CrsMatrix too. All Ifpack2
1115// preconditioners can and should do dynamic casts if they need a type
1116// more specific than RowMatrix. This keeps build time short and
1117// library and executable sizes small.
1118
1119#ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1120
1121#define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \
1122 template \
1123 class Ifpack2::BlockRelaxation< \
1124 Tpetra::RowMatrix<S, LO, GO, N>, \
1125 Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N> > >;
1126
1127#endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1128
1129#endif // IFPACK2_BLOCKRELAXATION_DEF_HPP
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:341
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:350
int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:398
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_BlockRelaxation_def.hpp:124
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:673
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_BlockRelaxation_def.hpp:429
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:1047
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:118
int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:390
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:422
int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:383
void initialize()
Initialize.
Definition: Ifpack2_BlockRelaxation_def.hpp:559
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:364
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:543
double getComputeTime() const
Returns the time spent in compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:414
void setParameters(const Teuchos::ParameterList &params)
Sets all the parameters for the preconditioner.
Definition: Ifpack2_BlockRelaxation_def.hpp:173
double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:406
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:446
std::string description() const
A one-line description of this object.
Definition: Ifpack2_BlockRelaxation_def.hpp:988
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition: Ifpack2_BlockRelaxation_def.hpp:327
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_BlockRelaxation_decl.hpp:108
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
Interface for creating and solving a set of local linear problems.
Definition: Ifpack2_Container_decl.hpp:113
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:77
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:73
void getValidParameters(Teuchos::ParameterList &params)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition: Ifpack2_Parameters.cpp:51
void getParameter(const Teuchos::ParameterList &params, 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