Ifpack2 Templated Preconditioning Package  Version 1.0
Ifpack2_BlockRelaxation_def.hpp
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Tempated 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 
53 namespace Ifpack2 {
54 
55 template<class MatrixType,class ContainerType>
57 setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
58 {
59  if (A.getRawPtr () != A_.getRawPtr ()) { // it's a different matrix
60  IsInitialized_ = false;
61  IsComputed_ = false;
62  Partitioner_ = Teuchos::null;
63  W_ = Teuchos::null;
64 
65  if (! A.is_null ()) {
66  IsParallel_ = (A->getRowMap ()->getComm ()->getSize () != 1);
67  Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
68  Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A);
69  hasBlockCrsMatrix_ = !A_bcrs.is_null();
70  }
71  if (! Container_.is_null ()) {
72  Container_->clearBlocks();
73  }
74  NumLocalBlocks_ = 0;
75 
76  A_ = A;
77  computeImporter();
78  }
79 }
80 
81 template<class MatrixType,class ContainerType>
83 BlockRelaxation (const Teuchos::RCP<const row_matrix_type>& A)
84 :
85  Time_ (Teuchos::rcp (new Teuchos::Time ("Ifpack2::BlockRelaxation"))),
86  Container_ (Teuchos::null),
87  Partitioner_ (Teuchos::null),
88  PartitionerType_ ("linear"),
89  NumSweeps_ (1),
90  NumLocalBlocks_(0),
91  containerType_ ("Dense"),
92  PrecType_ (Ifpack2::Details::JACOBI),
93  ZeroStartingSolution_ (true),
94  hasBlockCrsMatrix_ (false),
95  DoBackwardGS_ (false),
96  OverlapLevel_ (0),
97  DampingFactor_ (STS::one ()),
98  IsInitialized_ (false),
99  IsComputed_ (false),
100  NumInitialize_ (0),
101  NumCompute_ (0),
102  NumApply_ (0),
103  InitializeTime_ (0.0),
104  ComputeTime_ (0.0),
105  NumLocalRows_ (0),
106  NumGlobalRows_ (0),
107  NumGlobalNonzeros_ (0),
108  W_ (Teuchos::null),
109  Importer_ (Teuchos::null)
110 {
111  setMatrix(A);
112 }
113 
114 template<class MatrixType,class ContainerType>
117 {}
118 
119 template<class MatrixType,class ContainerType>
120 void
122 setParameters (const Teuchos::ParameterList& List)
123 {
124  // Note that the validation process does not change List.
125  Teuchos::ParameterList validparams;
126  Ifpack2::getValidParameters (validparams);
127  List.validateParameters (validparams);
128 
129  if (List.isParameter ("relaxation: container")) {
130  // If the container type isn't a string, this will throw, but it
131  // rightfully should.
132  containerType_ = List.get<std::string> ("relaxation: container");
133  }
134 
135  if (List.isParameter ("relaxation: type")) {
136  const std::string relaxType = List.get<std::string> ("relaxation: type");
137 
138  if (relaxType == "Jacobi") {
139  PrecType_ = Ifpack2::Details::JACOBI;
140  }
141  else if (relaxType == "Gauss-Seidel") {
142  PrecType_ = Ifpack2::Details::GS;
143  }
144  else if (relaxType == "Symmetric Gauss-Seidel") {
145  PrecType_ = Ifpack2::Details::SGS;
146  }
147  else {
148  TEUCHOS_TEST_FOR_EXCEPTION
149  (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
150  "setParameters: Invalid parameter value \"" << relaxType
151  << "\" for parameter \"relaxation: type\".");
152  }
153  }
154 
155  if (List.isParameter ("relaxation: sweeps")) {
156  NumSweeps_ = List.get<int> ("relaxation: sweeps");
157  }
158 
159  // Users may specify this as a floating-point literal. In that
160  // case, it may have type double, even if scalar_type is something
161  // else. We take care to try the various types that make sense.
162  if (List.isParameter ("relaxation: damping factor")) {
163  if (List.isType<double> ("relaxation: damping factor")) {
164  const double dampingFactor =
165  List.get<double> ("relaxation: damping factor");
166  DampingFactor_ = static_cast<scalar_type> (dampingFactor);
167  }
168  else if (List.isType<scalar_type> ("relaxation: damping factor")) {
169  DampingFactor_ = List.get<scalar_type> ("relaxation: damping factor");
170  }
171  else if (List.isType<magnitude_type> ("relaxation: damping factor")) {
172  const magnitude_type dampingFactor =
173  List.get<magnitude_type> ("relaxation: damping factor");
174  DampingFactor_ = static_cast<scalar_type> (dampingFactor);
175  }
176  else {
177  TEUCHOS_TEST_FOR_EXCEPTION
178  (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
179  "setParameters: Parameter \"relaxation: damping factor\" "
180  "has the wrong type.");
181  }
182  }
183 
184  if (List.isParameter ("relaxation: zero starting solution")) {
185  ZeroStartingSolution_ = List.get<bool> ("relaxation: zero starting solution");
186  }
187 
188  if (List.isParameter ("relaxation: backward mode")) {
189  DoBackwardGS_ = List.get<bool> ("relaxation: backward mode");
190  }
191 
192  if (List.isParameter ("partitioner: type")) {
193  PartitionerType_ = List.get<std::string> ("partitioner: type");
194  }
195 
196  // Users may specify this as an int literal, so we need to allow
197  // both int and local_ordinal_type (not necessarily same as int).
198  if (List.isParameter ("partitioner: local parts")) {
199  if (List.isType<local_ordinal_type> ("partitioner: local parts")) {
200  NumLocalBlocks_ = List.get<local_ordinal_type> ("partitioner: local parts");
201  }
202  else if (! std::is_same<int, local_ordinal_type>::value &&
203  List.isType<int> ("partitioner: local parts")) {
204  NumLocalBlocks_ = List.get<int> ("partitioner: local parts");
205  }
206  else {
207  TEUCHOS_TEST_FOR_EXCEPTION
208  (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
209  "setParameters: Parameter \"partitioner: local parts\" "
210  "has the wrong type.");
211  }
212  }
213 
214  if (List.isParameter ("partitioner: overlap level")) {
215  if (List.isType<int> ("partitioner: overlap level")) {
216  OverlapLevel_ = List.get<int> ("partitioner: overlap level");
217  }
218  else if (! std::is_same<int, local_ordinal_type>::value &&
219  List.isType<local_ordinal_type> ("partitioner: overlap level")) {
220  OverlapLevel_ = List.get<local_ordinal_type> ("partitioner: overlap level");
221  }
222  else {
223  TEUCHOS_TEST_FOR_EXCEPTION
224  (true, std::invalid_argument, "Ifpack2::BlockRelaxation::"
225  "setParameters: Parameter \"partitioner: overlap level\" "
226  "has the wrong type.");
227  }
228  }
229 
230  std::string defaultContainer = "Dense";
231  if(std::is_same<ContainerType, Container<MatrixType> >::value)
232  {
233  //Generic container template parameter, container type specified in List
234  Ifpack2::getParameter(List, "relaxation: container", defaultContainer);
235  }
236  // check parameters
237  if (PrecType_ != Ifpack2::Details::JACOBI) {
238  OverlapLevel_ = 0;
239  }
240  if (NumLocalBlocks_ < static_cast<local_ordinal_type> (0)) {
241  NumLocalBlocks_ = A_->getNodeNumRows() / (-NumLocalBlocks_);
242  }
243  // other checks are performed in Partitioner_
244 
245  // NTS: Sanity check to be removed at a later date when Backward mode is enabled
246  TEUCHOS_TEST_FOR_EXCEPTION(
247  DoBackwardGS_, std::runtime_error,
248  "Ifpack2::BlockRelaxation:setParameters: Setting the \"relaxation: "
249  "backward mode\" parameter to true is not yet supported.");
250 
251  // copy the list as each subblock's constructor will
252  // require it later
253  List_ = List;
254 }
255 
256 template<class MatrixType,class ContainerType>
257 Teuchos::RCP<const Teuchos::Comm<int> >
259 {
260  TEUCHOS_TEST_FOR_EXCEPTION
261  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::getComm: "
262  "The matrix is null. You must call setMatrix() with a nonnull matrix "
263  "before you may call this method.");
264  return A_->getComm ();
265 }
266 
267 template<class MatrixType,class ContainerType>
268 Teuchos::RCP<const Tpetra::RowMatrix<typename MatrixType::scalar_type,
269  typename MatrixType::local_ordinal_type,
270  typename MatrixType::global_ordinal_type,
271  typename MatrixType::node_type> >
273  return A_;
274 }
275 
276 template<class MatrixType,class ContainerType>
277 Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
278  typename MatrixType::global_ordinal_type,
279  typename MatrixType::node_type> >
281 getDomainMap () const
282 {
283  TEUCHOS_TEST_FOR_EXCEPTION
284  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
285  "getDomainMap: The matrix is null. You must call setMatrix() with a "
286  "nonnull matrix before you may call this method.");
287  return A_->getDomainMap ();
288 }
289 
290 template<class MatrixType,class ContainerType>
291 Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
292  typename MatrixType::global_ordinal_type,
293  typename MatrixType::node_type> >
295 getRangeMap () const
296 {
297  TEUCHOS_TEST_FOR_EXCEPTION
298  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
299  "getRangeMap: The matrix is null. You must call setMatrix() with a "
300  "nonnull matrix before you may call this method.");
301  return A_->getRangeMap ();
302 }
303 
304 template<class MatrixType,class ContainerType>
305 bool
307 hasTransposeApply () const {
308  return true;
309 }
310 
311 template<class MatrixType,class ContainerType>
312 int
315  return NumInitialize_;
316 }
317 
318 template<class MatrixType,class ContainerType>
319 int
322 {
323  return NumCompute_;
324 }
325 
326 template<class MatrixType,class ContainerType>
327 int
329 getNumApply () const
330 {
331  return NumApply_;
332 }
333 
334 template<class MatrixType,class ContainerType>
335 double
338 {
339  return InitializeTime_;
340 }
341 
342 template<class MatrixType,class ContainerType>
343 double
346 {
347  return ComputeTime_;
348 }
349 
350 template<class MatrixType,class ContainerType>
351 double
353 getApplyTime () const
354 {
355  return ApplyTime_;
356 }
357 
358 
359 template<class MatrixType,class ContainerType>
361  TEUCHOS_TEST_FOR_EXCEPTION(
362  A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::getNodeSmootherComplexity: "
363  "The input matrix A is null. Please call setMatrix() with a nonnull "
364  "input matrix, then call compute(), before calling this method.");
365  // Relaxation methods cost roughly one apply + one block-diagonal inverse per iteration
366  // NOTE: This approximates all blocks as dense, which may overstate the cost if you have a sparse (or banded) container.
367  size_t block_nnz = 0;
368  for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i)
369  block_nnz += Partitioner_->numRowsInPart(i) *Partitioner_->numRowsInPart(i);
370 
371  return block_nnz + A_->getNodeNumEntries();
372 }
373 
374 template<class MatrixType,class ContainerType>
375 void
377 apply (const Tpetra::MultiVector<typename MatrixType::scalar_type,
378  typename MatrixType::local_ordinal_type,
379  typename MatrixType::global_ordinal_type,
380  typename MatrixType::node_type>& X,
381  Tpetra::MultiVector<typename MatrixType::scalar_type,
382  typename MatrixType::local_ordinal_type,
383  typename MatrixType::global_ordinal_type,
384  typename MatrixType::node_type>& Y,
385  Teuchos::ETransp mode,
386  scalar_type alpha,
387  scalar_type beta) const
388 {
389  TEUCHOS_TEST_FOR_EXCEPTION
390  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::apply: "
391  "The matrix is null. You must call setMatrix() with a nonnull matrix, "
392  "then call initialize() and compute() (in that order), before you may "
393  "call this method.");
394  TEUCHOS_TEST_FOR_EXCEPTION(
395  ! isComputed (), std::runtime_error, "Ifpack2::BlockRelaxation::apply: "
396  "isComputed() must be true prior to calling apply.");
397  TEUCHOS_TEST_FOR_EXCEPTION(
398  X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
399  "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = "
400  << X.getNumVectors () << " != Y.getNumVectors() = "
401  << Y.getNumVectors () << ".");
402  TEUCHOS_TEST_FOR_EXCEPTION(
403  mode != Teuchos::NO_TRANS, std::logic_error, "Ifpack2::BlockRelaxation::"
404  "apply: This method currently only implements the case mode == Teuchos::"
405  "NO_TRANS. Transposed modes are not currently supported.");
406  TEUCHOS_TEST_FOR_EXCEPTION(
407  alpha != Teuchos::ScalarTraits<scalar_type>::one (), std::logic_error,
408  "Ifpack2::BlockRelaxation::apply: This method currently only implements "
409  "the case alpha == 1. You specified alpha = " << alpha << ".");
410  TEUCHOS_TEST_FOR_EXCEPTION(
411  beta != Teuchos::ScalarTraits<scalar_type>::zero (), std::logic_error,
412  "Ifpack2::BlockRelaxation::apply: This method currently only implements "
413  "the case beta == 0. You specified beta = " << beta << ".");
414 
415  Time_->start(true);
416 
417  // If X and Y are pointing to the same memory location,
418  // we need to create an auxiliary vector, Xcopy
419  Teuchos::RCP<const MV> X_copy;
420  {
421  auto X_lcl_host = X.template getLocalView<Kokkos::HostSpace> ();
422  auto Y_lcl_host = Y.template getLocalView<Kokkos::HostSpace> ();
423  if (X_lcl_host.ptr_on_device () == Y_lcl_host.ptr_on_device ()) {
424  X_copy = rcp (new MV (X, Teuchos::Copy));
425  } else {
426  X_copy = rcpFromRef (X);
427  }
428  }
429 
430  if (ZeroStartingSolution_) {
431  Y.putScalar (STS::zero ());
432  }
433 
434  // Flops are updated in each of the following.
435  switch (PrecType_) {
436  case Ifpack2::Details::JACOBI:
437  ApplyInverseJacobi(*X_copy,Y);
438  break;
439  case Ifpack2::Details::GS:
440  ApplyInverseGS(*X_copy,Y);
441  break;
442  case Ifpack2::Details::SGS:
443  ApplyInverseSGS(*X_copy,Y);
444  break;
445  default:
446  TEUCHOS_TEST_FOR_EXCEPTION
447  (true, std::logic_error, "Ifpack2::BlockRelaxation::apply: Invalid "
448  "PrecType_ enum value " << PrecType_ << ". Valid values are Ifpack2::"
449  "Details::JACOBI = " << Ifpack2::Details::JACOBI << ", Ifpack2::Details"
450  "::GS = " << Ifpack2::Details::GS << ", and Ifpack2::Details::SGS = "
451  << Ifpack2::Details::SGS << ". Please report this bug to the Ifpack2 "
452  "developers.");
453  }
454 
455  ++NumApply_;
456  Time_->stop();
457  ApplyTime_ += Time_->totalElapsedTime();
458 }
459 
460 template<class MatrixType,class ContainerType>
461 void
463 applyMat (const Tpetra::MultiVector<typename MatrixType::scalar_type,
464  typename MatrixType::local_ordinal_type,
465  typename MatrixType::global_ordinal_type,
466  typename MatrixType::node_type>& X,
467  Tpetra::MultiVector<typename MatrixType::scalar_type,
468  typename MatrixType::local_ordinal_type,
469  typename MatrixType::global_ordinal_type,
470  typename MatrixType::node_type>& Y,
471  Teuchos::ETransp mode) const
472 {
473  A_->apply (X, Y, mode);
474 }
475 
476 template<class MatrixType,class ContainerType>
477 void
480 {
481  using Teuchos::rcp;
482  typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>
483  row_graph_type;
484 
485  TEUCHOS_TEST_FOR_EXCEPTION
486  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::initialize: "
487  "The matrix is null. You must call setMatrix() with a nonnull matrix "
488  "before you may call this method.");
489 
490  // Check whether we have a BlockCrsMatrix
491  Teuchos::RCP<const block_crs_matrix_type> A_bcrs =
492  Teuchos::rcp_dynamic_cast<const block_crs_matrix_type> (A_);
493  hasBlockCrsMatrix_ = !A_bcrs.is_null();
494  if (A_bcrs.is_null ()) {
495  hasBlockCrsMatrix_ = false;
496  }
497  else {
498  hasBlockCrsMatrix_ = true;
499  }
500 
501  IsInitialized_ = false;
502  Time_->start (true);
503 
504  NumLocalRows_ = A_->getNodeNumRows ();
505  NumGlobalRows_ = A_->getGlobalNumRows ();
506  NumGlobalNonzeros_ = A_->getGlobalNumEntries ();
507 
508  // NTS: Will need to add support for Zoltan2 partitions later Also,
509  // will need a better way of handling the Graph typing issue.
510  // Especially with ordinal types w.r.t the container.
511  Partitioner_ = Teuchos::null;
512 
513  if (PartitionerType_ == "linear") {
514  Partitioner_ =
515  rcp (new Ifpack2::LinearPartitioner<row_graph_type> (A_->getGraph ()));
516  } else if (PartitionerType_ == "line") {
517  Partitioner_ =
519  } else if (PartitionerType_ == "user") {
520  Partitioner_ =
521  rcp (new Ifpack2::Details::UserPartitioner<row_graph_type> (A_->getGraph () ) );
522  } else {
523  // We should have checked for this in setParameters(), so it's a
524  // logic_error, not an invalid_argument or runtime_error.
525  TEUCHOS_TEST_FOR_EXCEPTION
526  (true, std::logic_error, "Ifpack2::BlockRelaxation::initialize: Unknown "
527  "partitioner type " << PartitionerType_ << ". Valid values are "
528  "\"linear\", \"line\", and \"user\".");
529  }
530 
531  // need to partition the graph of A
532  Partitioner_->setParameters (List_);
533  Partitioner_->compute ();
534 
535  // get actual number of partitions
536  NumLocalBlocks_ = Partitioner_->numLocalParts ();
537 
538  // Note: Unlike Ifpack, we'll punt on computing W_ until compute(), which is where
539  // we assume that the type of relaxation has been chosen.
540 
541  if (A_->getComm()->getSize() != 1) {
542  IsParallel_ = true;
543  } else {
544  IsParallel_ = false;
545  }
546 
547  ++NumInitialize_;
548  Time_->stop ();
549  InitializeTime_ += Time_->totalElapsedTime ();
550  IsInitialized_ = true;
551 }
552 
553 
554 template<class MatrixType,class ContainerType>
555 void
558 {
559  // typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type; // unused
560  // typedef Tpetra::Import<local_ordinal_type,global_ordinal_type, node_type> import_type; // unused
561  using Teuchos::Ptr;
562  using Teuchos::RCP;
563  using Teuchos::rcp;
564  using Teuchos::Array;
565  using Teuchos::ArrayView;
566 
567  Time_->start (true);
568 
569  // reset values
570  IsComputed_ = false;
571 
572  // Extract the submatrices
573  ExtractSubmatrices ();
574 
575  // Compute the weight vector if we're doing overlapped Jacobi (and
576  // only if we're doing overlapped Jacobi).
577  TEUCHOS_TEST_FOR_EXCEPTION
578  (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0, std::runtime_error,
579  "Ifpack2::BlockRelaxation::computeBlockCrs: "
580  "We do not support overlapped Jacobi yet for Tpetra::BlockCrsMatrix. Sorry!");
581 
582  ++NumCompute_;
583  Time_->stop ();
584  ComputeTime_ += Time_->totalElapsedTime();
585  IsComputed_ = true;
586 }
587 
588 
589 template<class MatrixType,class ContainerType>
590 void
593 {
594  using Teuchos::rcp;
595  typedef Tpetra::Vector<scalar_type,
597  // typedef Tpetra::Import<local_ordinal_type,
598  // global_ordinal_type, node_type> import_type; // unused
599 
600  TEUCHOS_TEST_FOR_EXCEPTION
601  (A_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::compute: "
602  "The matrix is null. You must call setMatrix() with a nonnull matrix, "
603  "then call initialize(), before you may call this method.");
604 
605  // We should have checked for this in setParameters(), so it's a
606  // logic_error, not an invalid_argument or runtime_error.
607  TEUCHOS_TEST_FOR_EXCEPTION
608  (NumSweeps_ < 0, std::logic_error, "Ifpack2::BlockRelaxation::compute: "
609  "NumSweeps_ = " << NumSweeps_ << " < 0.");
610 
611  if (! isInitialized ()) {
612  initialize ();
613  }
614 
615  if (hasBlockCrsMatrix_) {
616  computeBlockCrs ();
617  return;
618  }
619 
620  Time_->start (true);
621 
622  // reset values
623  IsComputed_ = false;
624 
625  // Extract the submatrices
626  ExtractSubmatrices ();
627 
628  // Compute the weight vector if we're doing overlapped Jacobi (and
629  // only if we're doing overlapped Jacobi).
630  if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) {
631  // weight of each vertex
632  W_ = rcp (new vector_type (A_->getRowMap ()));
633  W_->putScalar (STS::zero ());
634  Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0);
635 
636  for (local_ordinal_type i = 0 ; i < NumLocalBlocks_ ; ++i) {
637  for (size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) {
638  // FIXME (mfh 12 Sep 2014) Should this really be int?
639  // Perhaps it should be local_ordinal_type instead.
640  int LID = (*Partitioner_)(i,j);
641  w_ptr[LID]+= STS::one();
642  }
643  }
644  W_->reciprocal (*W_);
645  }
646  ++NumCompute_;
647  Time_->stop ();
648  ComputeTime_ += Time_->totalElapsedTime();
649  IsComputed_ = true;
650 }
651 
652 template<class MatrixType, class ContainerType>
653 void
656 {
657  TEUCHOS_TEST_FOR_EXCEPTION
658  (Partitioner_.is_null (), std::runtime_error, "Ifpack2::BlockRelaxation::"
659  "ExtractSubmatrices: Partitioner object is null.");
660 
661  NumLocalBlocks_ = Partitioner_->numLocalParts ();
662  std::string containerType = ContainerType::getName ();
663  if (containerType == "Generic") {
664  // ContainerType is Container<row_matrix_type>. Thus, we need to
665  // get the container name from the parameter list. We use "Dense"
666  // as the default container type.
667  containerType = containerType_;
668  }
669 
670  Teuchos::Array<Teuchos::Array<local_ordinal_type> > localRows(NumLocalBlocks_);
671  for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) {
672  const size_t numRows = Partitioner_->numRowsInPart (i);
673 
674  localRows[i].resize(numRows);
675  // Extract a list of the indices of each partitioner row.
676  for (size_t j = 0; j < numRows; ++j) {
677  localRows[i][j] = (*Partitioner_) (i,j);
678  }
679  }
680  Container_ = Teuchos::rcp_static_cast<Container<MatrixType>>
681  (Details::createContainer<row_matrix_type> (containerType, A_, localRows, Importer_,
682  OverlapLevel_, DampingFactor_));
683  Container_->setParameters(List_);
684  Container_->initialize();
685  Container_->compute(); //initialize + compute each block matrix
686 }
687 
688 template<class MatrixType,class ContainerType>
689 void
691 ApplyInverseJacobi (const MV& X, MV& Y) const
692 {
693  const size_t NumVectors = X.getNumVectors ();
694  typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace>();
695  typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace>();
696  MV AY (Y.getMap (), NumVectors);
697 
698  typename ContainerType::HostView AYView = AY.template getLocalView<Kokkos::HostSpace>();
699  // Initial matvec not needed
700  int starting_iteration = 0;
701  if (OverlapLevel_ > 0)
702  {
703  //Overlapping jacobi, with view of W_
704  typename ContainerType::HostView WView = W_->template getLocalView<Kokkos::HostSpace>();
705  if(ZeroStartingSolution_) {
706  Container_->DoOverlappingJacobi(XView, YView, WView, X.getStride());
707  starting_iteration = 1;
708  }
709  const scalar_type ONE = STS::one();
710  for(int j = starting_iteration; j < NumSweeps_; j++)
711  {
712  applyMat (Y, AY);
713  AY.update (ONE, X, -ONE);
714  Container_->DoOverlappingJacobi (AYView, YView, WView, X.getStride());
715  }
716  }
717  else
718  {
719  //Non-overlapping
720  if(ZeroStartingSolution_)
721  {
722  Container_->DoJacobi (XView, YView, X.getStride());
723  starting_iteration = 1;
724  }
725  const scalar_type ONE = STS::one();
726  for(int j = starting_iteration; j < NumSweeps_; j++)
727  {
728  applyMat (Y, AY);
729  AY.update (ONE, X, -ONE);
730  Container_->DoJacobi (AYView, YView, X.getStride());
731  }
732  }
733 }
734 
735 template<class MatrixType,class ContainerType>
736 void
738 ApplyInverseGS (const MV& X, MV& Y) const
739 {
740  using Teuchos::Ptr;
741  using Teuchos::ptr;
742  size_t numVecs = X.getNumVectors();
743  //Get view of X (is never modified in this function)
744  typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace>();
745  typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace>();
746  //Pre-import Y, if parallel
747  Ptr<MV> Y2;
748  bool deleteY2 = false;
749  if(IsParallel_)
750  {
751  Y2 = ptr(new MV(Importer_->getTargetMap(), numVecs));
752  deleteY2 = true;
753  }
754  else
755  Y2 = ptr(&Y);
756  if(IsParallel_)
757  {
758  for(int j = 0; j < NumSweeps_; ++j)
759  {
760  //do import once per sweep
761  Y2->doImport(Y, *Importer_, Tpetra::INSERT);
762  typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace>();
763  Container_->DoGaussSeidel(XView, YView, Y2View, X.getStride());
764  }
765  }
766  else
767  {
768  typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace>();
769  for(int j = 0; j < NumSweeps_; ++j)
770  {
771  Container_->DoGaussSeidel(XView, YView, Y2View, X.getStride());
772  }
773  }
774  if(deleteY2)
775  delete Y2.get();
776 }
777 
778 template<class MatrixType,class ContainerType>
779 void
781 ApplyInverseSGS (const MV& X, MV& Y) const
782 {
783  using Teuchos::Ptr;
784  using Teuchos::ptr;
785  //Get view of X (is never modified in this function)
786  typename ContainerType::HostView XView = X.template getLocalView<Kokkos::HostSpace>();
787  typename ContainerType::HostView YView = Y.template getLocalView<Kokkos::HostSpace>();
788  //Pre-import Y, if parallel
789  Ptr<MV> Y2;
790  bool deleteY2 = false;
791  if(IsParallel_)
792  {
793  Y2 = ptr(new MV(Importer_->getTargetMap(), X.getNumVectors()));
794  deleteY2 = true;
795  }
796  else
797  Y2 = ptr(&Y);
798  if(IsParallel_)
799  {
800  for(int j = 0; j < NumSweeps_; ++j)
801  {
802  //do import once per sweep
803  Y2->doImport(Y, *Importer_, Tpetra::INSERT);
804  typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace>();
805  Container_->DoSGS(XView, YView, Y2View, X.getStride());
806  }
807  }
808  else
809  {
810  typename ContainerType::HostView Y2View = Y2->template getLocalView<Kokkos::HostSpace>();
811  for(int j = 0; j < NumSweeps_; ++j)
812  {
813  Container_->DoSGS(XView, YView, Y2View, X.getStride());
814  }
815  }
816  if(deleteY2)
817  delete Y2.get();
818 }
819 
820 template<class MatrixType,class ContainerType>
822 {
823  using Teuchos::RCP;
824  using Teuchos::rcp;
825  using Teuchos::Ptr;
826  using Teuchos::ArrayView;
827  using Teuchos::Array;
828  using Teuchos::Comm;
829  using Teuchos::rcp_dynamic_cast;
830  if(IsParallel_)
831  {
832  if(hasBlockCrsMatrix_)
833  {
834  const RCP<const block_crs_matrix_type> bcrs =
835  rcp_dynamic_cast<const block_crs_matrix_type>(A_);
836  int bs_ = bcrs->getBlockSize();
837  RCP<const map_type> oldDomainMap = A_->getDomainMap();
838  RCP<const map_type> oldColMap = A_->getColMap();
839  // Because A is a block CRS matrix, import will not do what you think it does
840  // We have to construct the correct maps for it
841  global_size_t numGlobalElements = oldColMap->getGlobalNumElements() * bs_;
842  global_ordinal_type indexBase = oldColMap->getIndexBase();
843  RCP<const Comm<int>> comm = oldColMap->getComm();
844  ArrayView<const global_ordinal_type> oldColElements = oldColMap->getNodeElementList();
845  Array<global_ordinal_type> newColElements(bs_ * oldColElements.size());
846  for(int i = 0; i < oldColElements.size(); i++)
847  {
848  for(int j = 0; j < bs_; j++)
849  newColElements[i * bs_ + j] = oldColElements[i] * bs_ + j;
850  }
851  RCP<map_type> colMap(new map_type(numGlobalElements, newColElements, indexBase, comm));
852  // Create the importer
853  Importer_ = rcp(new import_type(oldDomainMap, colMap));
854  }
855  else if(!A_.is_null())
856  {
857  Importer_ = A_->getGraph()->getImporter();
858  if(Importer_.is_null())
859  Importer_ = rcp(new import_type(A_->getDomainMap(), A_->getColMap()));
860  }
861  }
862  //otherwise, Importer_ is not needed and is left NULL
863 }
864 
865 template<class MatrixType, class ContainerType>
866 std::string
868 description () const
869 {
870  std::ostringstream out;
871 
872  // Output is a valid YAML dictionary in flow style. If you don't
873  // like everything on a single line, you should call describe()
874  // instead.
875  out << "\"Ifpack2::BlockRelaxation\": {";
876  if (this->getObjectLabel () != "") {
877  out << "Label: \"" << this->getObjectLabel () << "\", ";
878  }
879  out << "Initialized: " << (isInitialized () ? "true" : "false") << ", ";
880  out << "Computed: " << (isComputed () ? "true" : "false") << ", ";
881  if (A_.is_null ()) {
882  out << "Matrix: null, ";
883  }
884  else {
885  out << "Matrix: not null"
886  << ", Global matrix dimensions: ["
887  << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "], ";
888  }
889 
890  // It's useful to print this instance's relaxation method. If you
891  // want more info than that, call describe() instead.
892  out << "\"relaxation: type\": ";
893  if (PrecType_ == Ifpack2::Details::JACOBI) {
894  out << "Block Jacobi";
895  } else if (PrecType_ == Ifpack2::Details::GS) {
896  out << "Block Gauss-Seidel";
897  } else if (PrecType_ == Ifpack2::Details::SGS) {
898  out << "Block Symmetric Gauss-Seidel";
899  } else {
900  out << "INVALID";
901  }
902 
903  out << ", overlap: " << OverlapLevel_;
904 
905  out << ", " << "sweeps: " << NumSweeps_ << ", "
906  << "damping factor: " << DampingFactor_ << ", ";
907 
908  std::string containerType = ContainerType::getName();
909  out << "block container: " << ((containerType == "Generic") ? containerType_ : containerType);
910 
911  out << "}";
912  return out.str();
913 }
914 
915 template<class MatrixType,class ContainerType>
916 void
918 describe (Teuchos::FancyOStream& out,
919  const Teuchos::EVerbosityLevel verbLevel) const
920 {
921  using std::endl;
922  using std::setw;
923  using Teuchos::TypeNameTraits;
924  using Teuchos::VERB_DEFAULT;
925  using Teuchos::VERB_NONE;
926  using Teuchos::VERB_LOW;
927  using Teuchos::VERB_MEDIUM;
928  using Teuchos::VERB_HIGH;
929  using Teuchos::VERB_EXTREME;
930 
931  Teuchos::EVerbosityLevel vl = verbLevel;
932  if (vl == VERB_DEFAULT) vl = VERB_LOW;
933  const int myImageID = A_->getComm()->getRank();
934 
935  // Convention requires that describe() begin with a tab.
936  Teuchos::OSTab tab (out);
937 
938  // none: print nothing
939  // low: print O(1) info from node 0
940  // medium:
941  // high:
942  // extreme:
943  if (vl != VERB_NONE && myImageID == 0) {
944  out << "Ifpack2::BlockRelaxation<"
945  << TypeNameTraits<MatrixType>::name () << ", "
946  << TypeNameTraits<ContainerType>::name () << " >:";
947  Teuchos::OSTab tab1 (out);
948 
949  if (this->getObjectLabel () != "") {
950  out << "label: \"" << this->getObjectLabel () << "\"" << endl;
951  }
952  out << "initialized: " << (isInitialized () ? "true" : "false") << endl
953  << "computed: " << (isComputed () ? "true" : "false") << endl;
954 
955  std::string precType;
956  if (PrecType_ == Ifpack2::Details::JACOBI) {
957  precType = "Block Jacobi";
958  } else if (PrecType_ == Ifpack2::Details::GS) {
959  precType = "Block Gauss-Seidel";
960  } else if (PrecType_ == Ifpack2::Details::SGS) {
961  precType = "Block Symmetric Gauss-Seidel";
962  }
963  out << "type: " << precType << endl
964  << "global number of rows: " << A_->getGlobalNumRows () << endl
965  << "global number of columns: " << A_->getGlobalNumCols () << endl
966  << "number of sweeps: " << NumSweeps_ << endl
967  << "damping factor: " << DampingFactor_ << endl
968  << "backwards mode: "
969  << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ? "true" : "false")
970  << endl
971  << "zero starting solution: "
972  << (ZeroStartingSolution_ ? "true" : "false") << endl;
973  }
974 }
975 
976 }//namespace Ifpack2
977 
978 
979 // Macro that does explicit template instantiation (ETI) for
980 // Ifpack2::BlockRelaxation. S, LO, GO, N correspond to the four
981 // template parameters of Ifpack2::Preconditioner and
982 // Tpetra::RowMatrix.
983 //
984 // We only instantiate for MatrixType = Tpetra::RowMatrix. There's no
985 // need to instantiate for Tpetra::CrsMatrix too. All Ifpack2
986 // preconditioners can and should do dynamic casts if they need a type
987 // more specific than RowMatrix. This keeps build time short and
988 // library and executable sizes small.
989 
990 #ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION
991 
992 #define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \
993  template \
994  class Ifpack2::BlockRelaxation< \
995  Tpetra::RowMatrix<S, LO, GO, N>, \
996  Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N> > >;
997 
998 #endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION
999 
1000 #endif // IFPACK2_BLOCKRELAXATION_DEF_HPP
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition: Ifpack2_BlockRelaxation_def.hpp:258
Ifpack2::BlockRelaxation class declaration.
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:353
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_BlockRelaxation_def.hpp:360
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:104
double getComputeTime() const
Returns the time spent in compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:345
int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:314
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:918
double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_BlockRelaxation_def.hpp:337
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner&#39;s constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:272
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:281
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
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:377
std::string description() const
A one-line description of this object.
Definition: Ifpack2_BlockRelaxation_def.hpp:868
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_BlockRelaxation_def.hpp:295
void setParameters(const Teuchos::ParameterList &params)
Sets all the parameters for the preconditioner.
Definition: Ifpack2_BlockRelaxation_def.hpp:122
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_BlockRelaxation_decl.hpp:216
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition: Ifpack2_BlockRelaxation_def.hpp:592
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition: Ifpack2_BlockRelaxation_decl.hpp:83
Ifpack2 implementation details.
Declaration of a user-defined partitioner in which the user can define a partition of the graph...
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:463
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_BlockRelaxation_decl.hpp:109
Tpetra::Import< local_ordinal_type, global_ordinal_type, node_type > import_type
Tpetra::Importer specialization for use with MatrixType and compatible MultiVectors.
Definition: Ifpack2_BlockRelaxation_decl.hpp:115
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_BlockRelaxation_def.hpp:57
virtual ~BlockRelaxation()
Destructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:116
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
void initialize()
Initialize.
Definition: Ifpack2_BlockRelaxation_def.hpp:479
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_BlockRelaxation_decl.hpp:226
MatrixType::node_type node_type
Node type of the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:106
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:101
int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_BlockRelaxation_def.hpp:321
Definition: Ifpack2_Container.hpp:761
int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_BlockRelaxation_def.hpp:329
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_BlockRelaxation_decl.hpp:98
Interface for creating and solving a local linear problem.
Definition: Ifpack2_Container.hpp:114
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
void getValidParameters(Teuchos::ParameterList &params)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition: Ifpack2_Parameters.cpp:50
Ifpack2::LinePartitioner: A class to define partitions into a set of lines.
Definition: Ifpack2_LinePartitioner_decl.hpp:77
A class to define linear partitions.
Definition: Ifpack2_LinearPartitioner_decl.hpp:60
BlockRelaxation(const Teuchos::RCP< const row_matrix_type > &Matrix)
Constructor.
Definition: Ifpack2_BlockRelaxation_def.hpp:83