Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
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 nonsymCombine_(0),
100 schwarzCombineMode_("ZERO"),
101 DampingFactor_ (STS::one ()),
102 IsInitialized_ (false),
103 IsComputed_ (false),
104 NumInitialize_ (0),
105 NumCompute_ (0),
106 TimerForApply_(true),
107 NumApply_ (0),
108 InitializeTime_ (0.0),
109 ComputeTime_ (0.0),
110 NumLocalRows_ (0),
111 NumGlobalRows_ (0),
112 NumGlobalNonzeros_ (0),
113 W_ (Teuchos::null),
114 Importer_ (Teuchos::null)
115{
116 setMatrix(A);
117}
118
119template<class MatrixType,class ContainerType>
123
124template<class MatrixType,class ContainerType>
125Teuchos::RCP<const Teuchos::ParameterList>
127getValidParameters () const
128{
129 Teuchos::RCP<Teuchos::ParameterList> validParams = Teuchos::parameterList ("Ifpack2::BlockRelaxation");
130
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); // mfh 24 Mar 2015: for backwards compatibility ONLY
139 validParams->set("schwarz: combine mode", "ZERO"); // use string mode for this
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"); // use string mode for this
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);
157
158 Teuchos::ParameterList dummyList;
159 validParams->set("Amesos2",dummyList);
160 validParams->sublist("Amesos2").disableRecursiveValidation();
161 validParams->set("Amesos2 solver name", "KLU2");
162
163 Teuchos::ArrayRCP<int> tmp;
164 validParams->set("partitioner: map", tmp);
165
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);
174
175 return validParams;
176}
177
178template<class MatrixType,class ContainerType>
179void
181setParameters (const Teuchos::ParameterList& pl)
182{
183 // CAG: Copied from Relaxation
184 // FIXME (aprokop 18 Oct 2013) Casting away const is bad here.
185 // but otherwise, we will get [unused] in pl
186 this->setParametersImpl(const_cast<Teuchos::ParameterList&>(pl));
187}
188
189template<class MatrixType,class ContainerType>
190void
192setParametersImpl (Teuchos::ParameterList& List)
193{
194 if (List.isType<double>("relaxation: damping factor")) {
195 // Make sure that ST=complex can run with a damping factor that is
196 // a double.
197 scalar_type df = List.get<double>("relaxation: damping factor");
198 List.remove("relaxation: damping factor");
199 List.set("relaxation: damping factor",df);
200 }
201
202 // Note that the validation process does not change List.
203 Teuchos::RCP<const Teuchos::ParameterList> validparams;
204 validparams = this->getValidParameters();
205 List.validateParameters (*validparams);
206
207 if (List.isParameter ("relaxation: container")) {
208 // If the container type isn't a string, this will throw, but it
209 // rightfully should.
210
211 // If its value does not match the currently registered Container types,
212 // the ContainerFactory will throw with an informative message.
213 containerType_ = List.get<std::string> ("relaxation: container");
214 // Intercept more human-readable aliases for the *TriDi containers
215 if(containerType_ == "Tridiagonal") {
216 containerType_ = "TriDi";
217 }
218 if(containerType_ == "Block Tridiagonal") {
219 containerType_ = "BlockTriDi";
220 }
221 }
222
223 if (List.isParameter ("relaxation: type")) {
224 const std::string relaxType = List.get<std::string> ("relaxation: type");
225
226 if (relaxType == "Jacobi") {
227 PrecType_ = Ifpack2::Details::JACOBI;
228 }
229 else if (relaxType == "MT Split Jacobi") {
230 PrecType_ = Ifpack2::Details::MTSPLITJACOBI;
231 }
232 else if (relaxType == "Gauss-Seidel") {
233 PrecType_ = Ifpack2::Details::GS;
234 }
235 else if (relaxType == "Symmetric Gauss-Seidel") {
236 PrecType_ = Ifpack2::Details::SGS;
237 }
238 else {
239 TEUCHOS_TEST_FOR_EXCEPTION
240 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
241 "setParameters: Invalid parameter value \"" << relaxType
242 << "\" for parameter \"relaxation: type\".");
243 }
244 }
245
246 if (List.isParameter ("relaxation: sweeps")) {
247 NumSweeps_ = List.get<int> ("relaxation: sweeps");
248 }
249
250 // Users may specify this as a floating-point literal. In that
251 // case, it may have type double, even if scalar_type is something
252 // else. We take care to try the various types that make sense.
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);
258 }
259 else if (List.isType<scalar_type> ("relaxation: damping factor")) {
260 DampingFactor_ = List.get<scalar_type> ("relaxation: damping factor");
261 }
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);
266 }
267 else {
268 TEUCHOS_TEST_FOR_EXCEPTION
269 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
270 "setParameters: Parameter \"relaxation: damping factor\" "
271 "has the wrong type.");
272 }
273 }
274
275 if (List.isParameter ("relaxation: zero starting solution")) {
276 ZeroStartingSolution_ = List.get<bool> ("relaxation: zero starting solution");
277 }
278
279 if (List.isParameter ("relaxation: backward mode")) {
280 DoBackwardGS_ = List.get<bool> ("relaxation: backward mode");
281 }
282
283 if (List.isParameter ("partitioner: type")) {
284 PartitionerType_ = List.get<std::string> ("partitioner: type");
285 }
286
287 // Users may specify this as an int literal, so we need to allow
288 // both int and local_ordinal_type (not necessarily same as int).
289 if (List.isParameter ("partitioner: local parts")) {
290 if (List.isType<local_ordinal_type> ("partitioner: local parts")) {
291 NumLocalBlocks_ = List.get<local_ordinal_type> ("partitioner: local parts");
292 }
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");
296 }
297 else {
298 TEUCHOS_TEST_FOR_EXCEPTION
299 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
300 "setParameters: Parameter \"partitioner: local parts\" "
301 "has the wrong type.");
302 }
303 }
304
305 if (List.isParameter ("partitioner: overlap level")) {
306 if (List.isType<int> ("partitioner: overlap level")) {
307 OverlapLevel_ = List.get<int> ("partitioner: overlap level");
308 }
309 else if (! std::is_same<int, local_ordinal_type>::value &&
310 List.isType<local_ordinal_type> ("partitioner: overlap level")) {
311 OverlapLevel_ = List.get<local_ordinal_type> ("partitioner: overlap level");
312 }
313 else {
314 TEUCHOS_TEST_FOR_EXCEPTION
315 (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
316 "setParameters: Parameter \"partitioner: overlap level\" "
317 "has the wrong type.");
318 }
319 }
320 // when using global ID parts, assume that some blocks overlap even if
321 // user did not explicitly set the overlap level in the input file.
322 if ( ( List.isParameter("partitioner: global ID parts")) && (OverlapLevel_ < 1)) OverlapLevel_ = 1;
323
324 if (List.isParameter ("partitioner: nonsymmetric overlap combine"))
325 nonsymCombine_ = List.get<bool> ("partitioner: nonsymmetric overlap combine");
326
327 if (List.isParameter ("partitioner: combine mode"))
328 schwarzCombineMode_ = List.get<std::string> ("partitioner: combine mode");
329
330 std::string defaultContainer = "TriDi";
331 if(std::is_same<ContainerType, Container<MatrixType> >::value)
332 {
333 //Generic container template parameter, container type specified in List
334 Ifpack2::getParameter(List, "relaxation: container", defaultContainer);
335 }
336 // check parameters
337 if (PrecType_ != Ifpack2::Details::JACOBI) {
338 OverlapLevel_ = 0;
339 }
340 if (NumLocalBlocks_ < static_cast<local_ordinal_type> (0)) {
341 NumLocalBlocks_ = A_->getLocalNumRows() / (-NumLocalBlocks_);
342 }
343
344 decouple_ = false;
345 if(List.isParameter("block relaxation: decouple dofs"))
346 decouple_ = List.get<bool>("block relaxation: decouple dofs");
347
348 // other checks are performed in Partitioner_
349
350 // NTS: Sanity check to be removed at a later date when Backward mode is enabled
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.");
355
356 if(List.isParameter("timer for apply"))
357 TimerForApply_ = List.get<bool>("timer for apply");
358
359 // copy the list as each subblock's constructor will
360 // require it later
361 List_ = List;
362}
363
364template<class MatrixType,class ContainerType>
365Teuchos::RCP<const Teuchos::Comm<int> >
367{
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 ();
373}
374
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> >
383
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> >
389getDomainMap () const
390{
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 ();
396}
397
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> >
403getRangeMap () const
404{
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 ();
410}
411
412template<class MatrixType,class ContainerType>
413bool
415hasTransposeApply () const {
416 return true;
417}
418
419template<class MatrixType,class ContainerType>
420int
422getNumInitialize () const {
423 return NumInitialize_;
424}
425
426template<class MatrixType,class ContainerType>
427int
429getNumCompute () const
430{
431 return NumCompute_;
432}
433
434template<class MatrixType,class ContainerType>
435int
437getNumApply () const
438{
439 return NumApply_;
440}
441
442template<class MatrixType,class ContainerType>
443double
445getInitializeTime () const
446{
447 return InitializeTime_;
448}
449
450template<class MatrixType,class ContainerType>
451double
453getComputeTime () const
454{
455 return ComputeTime_;
456}
457
458template<class MatrixType,class ContainerType>
459double
461getApplyTime () const
462{
463 return ApplyTime_;
464}
465
466
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.");
473 // Relaxation methods cost roughly one apply + one block-diagonal inverse per iteration
474 // NOTE: This approximates all blocks as dense, which may overstate the cost if you have a sparse (or banded) container.
475 size_t block_nnz = 0;
476 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i)
477 block_nnz += Partitioner_->numRowsInPart(i) *Partitioner_->numRowsInPart(i);
478
479 return block_nnz + A_->getLocalNumEntries();
480}
481
482template<class MatrixType,class ContainerType>
483void
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,
494 scalar_type alpha,
495 scalar_type beta) const
496{
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 << ".");
522
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);
529 }
530 }
531
532 Teuchos::Time time = Teuchos::Time(timerName);
533 double startTime = time.wallTime();
534
535 {
536 Teuchos::RCP<Teuchos::TimeMonitor> timeMon;
537 if (TimerForApply_)
538 timeMon = Teuchos::rcp(new Teuchos::TimeMonitor(*timer));
539
540 // If X and Y are pointing to the same memory location,
541 // we need to create an auxiliary vector, Xcopy
542 Teuchos::RCP<const MV> X_copy;
543 {
544 if (X.aliases(Y)) {
545 X_copy = rcp (new MV (X, Teuchos::Copy));
546 } else {
547 X_copy = rcpFromRef (X);
548 }
549 }
550
551 if (ZeroStartingSolution_) {
552 Y.putScalar (STS::zero ());
553 }
554
555 // Flops are updated in each of the following.
556 switch (PrecType_) {
557 case Ifpack2::Details::JACOBI:
558 ApplyInverseJacobi(*X_copy,Y);
559 break;
560 case Ifpack2::Details::GS:
561 ApplyInverseGS(*X_copy,Y);
562 break;
563 case Ifpack2::Details::SGS:
564 ApplyInverseSGS(*X_copy,Y);
565 break;
566 case Ifpack2::Details::MTSPLITJACOBI:
567 //note: for this method, the container is always BlockTriDi
568 Container_->applyInverseJacobi(*X_copy, Y, DampingFactor_, ZeroStartingSolution_, NumSweeps_);
569 break;
570 default:
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 "
577 "developers.");
578 }
579 }
580
581 ApplyTime_ += (time.wallTime() - startTime);
582 ++NumApply_;
583}
584
585template<class MatrixType,class ContainerType>
586void
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
597{
598 A_->apply (X, Y, mode);
599}
600
601template<class MatrixType,class ContainerType>
602void
604initialize ()
605{
606 using Teuchos::rcp;
607 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
608 row_graph_type;
609
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.");
614
615 Teuchos::RCP<Teuchos::Time> timer =
616 Teuchos::TimeMonitor::getNewCounter ("Ifpack2::BlockRelaxation::initialize");
617 double startTime = timer->wallTime();
618
619 { // Timing of initialize starts here
620 Teuchos::TimeMonitor timeMon (*timer);
621 IsInitialized_ = false;
622
623 // Check whether we have a BlockCrsMatrix
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;
629 }
630 else {
631 hasBlockCrsMatrix_ = true;
632 }
633
634 NumLocalRows_ = A_->getLocalNumRows ();
635 NumGlobalRows_ = A_->getGlobalNumRows ();
636 NumGlobalNonzeros_ = A_->getGlobalNumEntries ();
637
638 // NTS: Will need to add support for Zoltan2 partitions later Also,
639 // will need a better way of handling the Graph typing issue.
640 // Especially with ordinal types w.r.t the container.
641 Partitioner_ = Teuchos::null;
642
643 if (PartitionerType_ == "linear") {
644 Partitioner_ =
645 rcp (new Ifpack2::LinearPartitioner<row_graph_type> (A_->getGraph ()));
646 } else if (PartitionerType_ == "line") {
647 Partitioner_ =
649 } else if (PartitionerType_ == "user") {
650 Partitioner_ =
651 rcp (new Ifpack2::Details::UserPartitioner<row_graph_type> (A_->getGraph () ) );
652 } else {
653 // We should have checked for this in setParameters(), so it's a
654 // logic_error, not an invalid_argument or runtime_error.
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\".");
659 }
660
661 // need to partition the graph of A
662 Partitioner_->setParameters (List_);
663 Partitioner_->compute ();
664
665 // get actual number of partitions
666 NumLocalBlocks_ = Partitioner_->numLocalParts ();
667
668 // Note: Unlike Ifpack, we'll punt on computing W_ until compute(), which is where
669 // we assume that the type of relaxation has been chosen.
670
671 if (A_->getComm()->getSize() != 1) {
672 IsParallel_ = true;
673 } else {
674 IsParallel_ = false;
675 }
676
677 // We should have checked for this in setParameters(), so it's a
678 // logic_error, not an invalid_argument or runtime_error.
679 TEUCHOS_TEST_FOR_EXCEPTION
680 (NumSweeps_ < 0, std::logic_error, "Ifpack2::BlockRelaxation::initialize: "
681 "NumSweeps_ = " << NumSweeps_ << " < 0.");
682
683 // Extract the submatrices
684 ExtractSubmatricesStructure ();
685
686 // Compute the weight vector if we're doing overlapped Jacobi (and
687 // only if we're doing overlapped Jacobi).
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!");
693
694 // weight of each vertex
695 W_ = rcp (new vector_type (A_->getRowMap ()));
696 W_->putScalar (STS::zero ());
697 {
698 Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0);
699
700 for (local_ordinal_type i = 0 ; i < NumLocalBlocks_ ; ++i) {
701 for (size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) {
702 local_ordinal_type LID = (*Partitioner_)(i,j);
703 w_ptr[LID] += STS::one();
704 }
705 }
706 }
707 // communicate to sum together W_[k]'s (# of blocks/patches) that update
708 // kth dof) and have this information in overlapped/extended vector.
709 // only needed when Schwarz combine mode is ADD as opposed to ZERO (which is RAS)
710
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);
724 }
725
726 }
727 W_->reciprocal (*W_);
728 }
729 } // timing of initialize stops here
730
731 InitializeTime_ += (timer->wallTime() - startTime);
732 ++NumInitialize_;
733 IsInitialized_ = true;
734}
735
736
737template<class MatrixType,class ContainerType>
738void
740compute ()
741{
742 using Teuchos::rcp;
743
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.");
748
749 if (! isInitialized ()) {
750 initialize ();
751 }
752
753 Teuchos::RCP<Teuchos::Time> timer =
754 Teuchos::TimeMonitor::getNewCounter ("Ifpack2::BlockRelaxation::compute");
755
756 double startTime = timer->wallTime();
757
758 {
759 Teuchos::TimeMonitor timeMon (*timer);
760
761 // reset values
762 IsComputed_ = false;
763
764 Container_->compute(); // compute each block matrix
765 }
766
767 ComputeTime_ += (timer->wallTime() - startTime);
768 ++NumCompute_;
769 IsComputed_ = true;
770}
771
772template<class MatrixType, class ContainerType>
773void
776{
777 TEUCHOS_TEST_FOR_EXCEPTION
778 (Partitioner_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
779 "ExtractSubmatricesStructure: Partitioner object is null.");
780
781 std::string containerType = ContainerType::getName ();
782 if (containerType == "Generic") {
783 // ContainerType is Container<row_matrix_type>. Thus, we need to
784 // get the container name from the parameter list. We use "TriDi"
785 // as the default container type.
786 containerType = containerType_;
787 }
788 //Whether the Container will define blocks (partitions)
789 //in terms of individual DOFs, and not by nodes (blocks).
790 bool pointIndexed = decouple_ && hasBlockCrsMatrix_;
791 Teuchos::Array<Teuchos::Array<local_ordinal_type> > blockRows;
792 if(decouple_)
793 {
794 //dofs [per node] is how many blocks each partition will be split into
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_)
800 {
801 for(local_ordinal_type i = 0; i < NumLocalBlocks_; i++)
802 {
803 size_t blockSize = Partitioner_->numRowsInPart(i);
804 //block i will be split into j different blocks,
805 //each corresponding to a different dof
806 for(local_ordinal_type j = 0; j < dofs; j++)
807 {
808 local_ordinal_type blockIndex = i * dofs + j;
809 blockRows[blockIndex].resize(blockSize);
810 for(size_t k = 0; k < blockSize; k++)
811 {
812 //here, the row and dof are combined to get the point index
813 //(what the row would be if A were a CrsMatrix)
814 blockRows[blockIndex][k] = (*Partitioner_)(i, k) * dofs + j;
815 }
816 }
817 }
818 }
819 else
820 {
821 //Rows in each partition are distributed round-robin to the blocks -
822 //that's how MueLu stores DOFs in a non-block matrix
823 for(local_ordinal_type i = 0; i < NumLocalBlocks_; i++)
824 {
825 //#ifdef HAVE_IFPACK2_DEBUG
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;
829 //#endif
830 //block i will be split into j different blocks,
831 //each corresponding to a different dof
832 for(local_ordinal_type j = 0; j < dofs; j++)
833 {
834 local_ordinal_type blockIndex = i * dofs + j;
835 blockRows[blockIndex].resize(blockSize);
836 for(size_t k = 0; k < blockSize; k++)
837 {
838 blockRows[blockIndex][k] = (*Partitioner_)(i, k * dofs + j);
839 }
840 }
841 }
842 }
843 }
844 else
845 {
846 //No decoupling - just get the rows directly from Partitioner
847 blockRows.resize(NumLocalBlocks_);
848 for(local_ordinal_type i = 0; i < NumLocalBlocks_; ++i)
849 {
850 const size_t numRows = Partitioner_->numRowsInPart (i);
851 blockRows[i].resize(numRows);
852 // Extract a list of the indices of each partitioner row.
853 for (size_t j = 0; j < numRows; ++j)
854 {
855 blockRows[i][j] = (*Partitioner_) (i,j);
856 }
857 }
858 }
859 //right before constructing the
860 Container_ = ContainerFactory<MatrixType>::build(containerType, A_, blockRows, Importer_, pointIndexed);
861 Container_->setParameters(List_);
862 Container_->initialize();
863}
864
865template<class MatrixType,class ContainerType>
866void
868ApplyInverseJacobi (const MV& X, MV& Y) const
869{
870 const size_t NumVectors = X.getNumVectors ();
871
872 MV AY (Y.getMap (), NumVectors);
873
874 // Initial matvec not needed
875 int starting_iteration = 0;
876 if (OverlapLevel_ > 0)
877 {
878 //Overlapping jacobi, with view of W_
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;
885 }
886 const scalar_type ONE = STS::one();
887 for(int j = starting_iteration; j < NumSweeps_; j++)
888 {
889 applyMat (Y, AY);
890 AY.update (ONE, X, -ONE);
891 {
892 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
893 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
894 Container_->DoOverlappingJacobi (AYView, YView, WView, DampingFactor_, nonsymCombine_);
895 }
896 }
897 }
898 else
899 {
900 //Non-overlapping
901 if(ZeroStartingSolution_)
902 {
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;
907 }
908 const scalar_type ONE = STS::one();
909 for(int j = starting_iteration; j < NumSweeps_; j++)
910 {
911 applyMat (Y, AY);
912 AY.update (ONE, X, -ONE);
913 {
914 auto AYView = AY.getLocalViewHost (Tpetra::Access::ReadOnly);
915 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
916 Container_->DoJacobi (AYView, YView, DampingFactor_);
917 }
918 }
919 }
920}
921
922template<class MatrixType,class ContainerType>
923void
925ApplyInverseGS (const MV& X, MV& Y) const
926{
927 using Teuchos::Ptr;
928 using Teuchos::ptr;
929 size_t numVecs = X.getNumVectors();
930 //Get view of X (is never modified in this function)
931 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
932 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
933 //Pre-import Y, if parallel
934 Ptr<MV> Y2;
935 bool deleteY2 = false;
936 if(IsParallel_)
937 {
938 Y2 = ptr(new MV(Importer_->getTargetMap(), numVecs));
939 deleteY2 = true;
940 }
941 else
942 Y2 = ptr(&Y);
943 if(IsParallel_)
944 {
945 for(int j = 0; j < NumSweeps_; ++j)
946 {
947 //do import once per sweep
948 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
949 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
950 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
951 }
952 }
953 else
954 {
955 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
956 for(int j = 0; j < NumSweeps_; ++j)
957 {
958 Container_->DoGaussSeidel(XView, YView, Y2View, DampingFactor_);
959 }
960 }
961 if(deleteY2)
962 delete Y2.get();
963}
964
965template<class MatrixType,class ContainerType>
966void
968ApplyInverseSGS (const MV& X, MV& Y) const
969{
970 using Teuchos::Ptr;
971 using Teuchos::ptr;
972 //Get view of X (is never modified in this function)
973 auto XView = X.getLocalViewHost (Tpetra::Access::ReadOnly);
974 auto YView = Y.getLocalViewHost (Tpetra::Access::ReadWrite);
975 //Pre-import Y, if parallel
976 Ptr<MV> Y2;
977 bool deleteY2 = false;
978 if(IsParallel_)
979 {
980 Y2 = ptr(new MV(Importer_->getTargetMap(), X.getNumVectors()));
981 deleteY2 = true;
982 }
983 else
984 Y2 = ptr(&Y);
985 if(IsParallel_)
986 {
987 for(int j = 0; j < NumSweeps_; ++j)
988 {
989 //do import once per sweep
990 Y2->doImport(Y, *Importer_, Tpetra::INSERT);
991 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
992 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
993 }
994 }
995 else
996 {
997 auto Y2View = Y2->getLocalViewHost (Tpetra::Access::ReadWrite);
998 for(int j = 0; j < NumSweeps_; ++j)
999 {
1000 Container_->DoSGS(XView, YView, Y2View, DampingFactor_);
1001 }
1002 }
1003 if(deleteY2)
1004 delete Y2.get();
1005}
1006
1007template<class MatrixType,class ContainerType>
1008void BlockRelaxation<MatrixType,ContainerType>::computeImporter () const
1009{
1010 using Teuchos::RCP;
1011 using Teuchos::rcp;
1012 using Teuchos::Ptr;
1013 using Teuchos::ArrayView;
1014 using Teuchos::Array;
1015 using Teuchos::Comm;
1016 using Teuchos::rcp_dynamic_cast;
1017 if(IsParallel_)
1018 {
1019 if(hasBlockCrsMatrix_)
1020 {
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();
1026 // Because A is a block CRS matrix, import will not do what you think it does
1027 // We have to construct the correct maps for it
1028 global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
1029 global_ordinal_type indexBase = oldColMap->getIndexBase();
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++)
1034 {
1035 for(int j = 0; j < bs_; j++)
1036 newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
1037 }
1038 RCP<map_type> colMap(new map_type(numGlobalElements, newColElements, indexBase, comm));
1039 // Create the importer
1040 Importer_ = rcp(new import_type(oldDomainMap, colMap));
1041 }
1042 else if(!A_.is_null())
1043 {
1044 Importer_ = A_->getGraph()->getImporter();
1045 if(Importer_.is_null())
1046 Importer_ = rcp(new import_type(A_->getDomainMap(), A_->getColMap()));
1047 }
1048 }
1049 //otherwise, Importer_ is not needed and is left NULL
1050}
1051
1052template<class MatrixType, class ContainerType>
1053std::string
1055description () const
1056{
1057 std::ostringstream out;
1058
1059 // Output is a valid YAML dictionary in flow style. If you don't
1060 // like everything on a single line, you should call describe()
1061 // instead.
1062 out << "\"Ifpack2::BlockRelaxation\": {";
1063 if (this->getObjectLabel () != "") {
1064 out << "Label: \"" << this->getObjectLabel () << "\", ";
1065 }
1066 // out << "Initialized: " << (isInitialized () ? "true" : "false") << ", ";
1067 // out << "Computed: " << (isComputed () ? "true" : "false") << ", ";
1068 if (A_.is_null ()) {
1069 out << "Matrix: null, ";
1070 }
1071 else {
1072 // out << "Matrix: not null"
1073 // << ", Global matrix dimensions: ["
1074 out << "Global matrix dimensions: ["
1075 << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "], ";
1076 }
1077
1078 // It's useful to print this instance's relaxation method. If you
1079 // want more info than that, call describe() instead.
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";
1089 } else {
1090 out << "INVALID";
1091 }
1092
1093 // BlockCrs if we have that
1094 if(hasBlockCrsMatrix_)
1095 out<<", BlockCrs";
1096
1097 // Print the approximate # rows per part
1098 int approx_rows_per_part = A_->getLocalNumRows()/Partitioner_->numLocalParts();
1099 out <<", blocksize: "<<approx_rows_per_part;
1100
1101 out << ", overlap: " << OverlapLevel_;
1102
1103 out << ", " << "sweeps: " << NumSweeps_ << ", "
1104 << "damping factor: " << DampingFactor_ << ", ";
1105
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");
1110
1111
1112 out << "}";
1113 return out.str();
1114}
1115
1116template<class MatrixType,class ContainerType>
1117void
1119describe (Teuchos::FancyOStream& out,
1120 const Teuchos::EVerbosityLevel verbLevel) const
1121{
1122 using std::endl;
1123 using std::setw;
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;
1131
1132 Teuchos::EVerbosityLevel vl = verbLevel;
1133 if (vl == VERB_DEFAULT) vl = VERB_LOW;
1134 const int myImageID = A_->getComm()->getRank();
1135
1136 // Convention requires that describe() begin with a tab.
1137 Teuchos::OSTab tab (out);
1138
1139 // none: print nothing
1140 // low: print O(1) info from node 0
1141 // medium:
1142 // high:
1143 // extreme:
1144 if (vl != VERB_NONE && myImageID == 0) {
1145 out << "Ifpack2::BlockRelaxation<"
1146 << TypeNameTraits<MatrixType>::name () << ", "
1147 << TypeNameTraits<ContainerType>::name () << " >:";
1148 Teuchos::OSTab tab1 (out);
1149
1150 if (this->getObjectLabel () != "") {
1151 out << "label: \"" << this->getObjectLabel () << "\"" << endl;
1152 }
1153 out << "initialized: " << (isInitialized () ? "true" : "false") << endl
1154 << "computed: " << (isComputed () ? "true" : "false") << endl;
1155
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";
1163 }
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")
1172 << endl
1173 << "zero starting solution: "
1174 << (ZeroStartingSolution_ ? "true" : "false") << endl;
1175 }
1176}
1177
1178}//namespace Ifpack2
1179
1180
1181// Macro that does explicit template instantiation (ETI) for
1182// Ifpack2::BlockRelaxation. S, LO, GO, N correspond to the four
1183// template parameters of Ifpack2::Preconditioner and
1184// Tpetra::RowMatrix.
1185//
1186// We only instantiate for MatrixType = Tpetra::RowMatrix. There's no
1187// need to instantiate for Tpetra::CrsMatrix too. All Ifpack2
1188// preconditioners can and should do dynamic casts if they need a type
1189// more specific than RowMatrix. This keeps build time short and
1190// library and executable sizes small.
1191
1192#ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1193
1194#define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \
1195 template \
1196 class Ifpack2::BlockRelaxation< \
1197 Tpetra::RowMatrix<S, LO, GO, N>, \
1198 Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N> > >;
1199
1200#endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION
1201
1202#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: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 &params)
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 &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