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