Ifpack2 Templated Preconditioning Package  Version 1.0
Ifpack2_DenseContainer_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_DENSECONTAINER_DEF_HPP
44 #define IFPACK2_DENSECONTAINER_DEF_HPP
45 
46 #include "Tpetra_CrsMatrix.hpp"
47 #include "Teuchos_LAPACK.hpp"
48 
49 #ifdef HAVE_MPI
50 # include <mpi.h>
51 # include "Teuchos_DefaultMpiComm.hpp"
52 #else
53 # include "Teuchos_DefaultSerialComm.hpp"
54 #endif // HAVE_MPI
55 
56 
57 namespace Ifpack2 {
58 
59 template<class MatrixType, class LocalScalarType>
61 DenseContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
62  const Teuchos::ArrayView<const local_ordinal_type>& localRows) :
63  Container<MatrixType> (matrix, localRows),
64  numRows_ (localRows.size ()),
65  diagBlock_ (numRows_, numRows_),
66  ipiv_ (numRows_, 0)
67 {
68  using Teuchos::Array;
69  using Teuchos::ArrayView;
70  using Teuchos::RCP;
71  using Teuchos::rcp;
72  using Teuchos::toString;
73  typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
74  typedef typename ArrayView<const local_ordinal_type>::size_type size_type;
75  TEUCHOS_TEST_FOR_EXCEPTION(
76  ! matrix->hasColMap (), std::invalid_argument, "Ifpack2::DenseContainer: "
77  "The constructor's input matrix must have a column Map.");
78 
79  // Check whether the input set of local row indices is correct.
80  const map_type& rowMap = * (matrix->getRowMap ());
81  const size_type numRows = localRows.size ();
82  bool rowIndicesValid = true;
83  Array<local_ordinal_type> invalidLocalRowIndices;
84  for (size_type i = 0; i < numRows; ++i) {
85  if (! rowMap.isNodeLocalElement (localRows[i])) {
86  rowIndicesValid = false;
87  invalidLocalRowIndices.push_back (localRows[i]);
88  break;
89  }
90  }
91  TEUCHOS_TEST_FOR_EXCEPTION(
92  ! rowIndicesValid, std::invalid_argument, "Ifpack2::DenseContainer: "
93  "On process " << rowMap.getComm ()->getRank () << " of "
94  << rowMap.getComm ()->getSize () << ", in the given set of local row "
95  "indices localRows = " << toString (localRows) << ", the following "
96  "entries are not valid local row indices on the calling process: "
97  << toString (invalidLocalRowIndices) << ".");
98 
99 #ifdef HAVE_MPI
100  RCP<const Teuchos::Comm<int> > localComm =
101  rcp (new Teuchos::MpiComm<int> (MPI_COMM_SELF));
102 #else
103  RCP<const Teuchos::Comm<int> > localComm =
104  rcp (new Teuchos::SerialComm<int> ());
105 #endif // HAVE_MPI
106 
107  // FIXME (mfh 25 Aug 2013) What if the matrix's row Map has a
108  // different index base than zero?
109  const global_ordinal_type indexBase = 0;
110  localMap_ = rcp (new map_type (numRows_, indexBase, localComm));
111 }
112 
113 template<class MatrixType, class LocalScalarType>
115 DenseContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
116  const Teuchos::ArrayView<const local_ordinal_type>& localRows) :
117  Container<MatrixType> (matrix, localRows)
118 {
119  TEUCHOS_TEST_FOR_EXCEPTION
120  (true, std::logic_error, "Ifpack2::DenseContainer: Not implemented for "
121  "LocalScalarType = " << Teuchos::TypeNameTraits<LocalScalarType>::name ()
122  << ".");
123 }
124 
125 template<class MatrixType, class LocalScalarType>
127 {}
128 
129 template<class MatrixType, class LocalScalarType>
131 {}
132 
133 
134 template<class MatrixType, class LocalScalarType>
136 setParameters (const Teuchos::ParameterList& /* List */)
137 {
138  // the solver doesn't currently take any parameters
139 }
140 
141 template<class MatrixType, class LocalScalarType>
143 setParameters (const Teuchos::ParameterList& /* List */)
144 {
145  // the solver doesn't currently take any parameters
146 }
147 
148 
149 template<class MatrixType, class LocalScalarType>
150 void
153 {
154  using Teuchos::null;
155  using Teuchos::rcp;
156 
157  // We assume that if you called this method, you intend to recompute
158  // everything.
159  IsInitialized_ = false;
160  IsComputed_ = false;
161 
162  // Fill the diagonal block and LU permutation array with zeros.
163  diagBlock_.putScalar (Teuchos::ScalarTraits<local_scalar_type>::zero ());
164  std::fill (ipiv_.begin (), ipiv_.end (), 0);
165 
166  IsInitialized_ = true;
167 }
168 
169 template<class MatrixType, class LocalScalarType>
170 void
173 {}
174 
175 template<class MatrixType, class LocalScalarType>
176 void
179 {
180  TEUCHOS_TEST_FOR_EXCEPTION(
181  static_cast<size_t> (ipiv_.size ()) != numRows_, std::logic_error,
182  "Ifpack2::DenseContainer::compute: ipiv_ array has the wrong size. "
183  "Please report this bug to the Ifpack2 developers.");
184 
185  IsComputed_ = false;
186  if (! this->isInitialized ()) {
187  this->initialize ();
188  }
189 
190  // Extract the submatrix.
191  extract (this->getMatrix ()); // extract the submatrix
192  factor (); // factor the submatrix
193 
194  IsComputed_ = true;
195 }
196 
197 template<class MatrixType, class LocalScalarType>
198 void
201 {}
202 
203 template<class MatrixType, class LocalScalarType>
204 void
206 factor ()
207 {
208  Teuchos::LAPACK<int, local_scalar_type> lapack;
209  int INFO = 0;
210  lapack.GETRF (diagBlock_.numRows (), diagBlock_.numCols (),
211  diagBlock_.values (), diagBlock_.stride (),
212  ipiv_.getRawPtr (), &INFO);
213  // INFO < 0 is a bug.
214  TEUCHOS_TEST_FOR_EXCEPTION(
215  INFO < 0, std::logic_error, "Ifpack2::DenseContainer::factor: "
216  "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
217  "incorrectly. INFO = " << INFO << " < 0. "
218  "Please report this bug to the Ifpack2 developers.");
219  // INFO > 0 means the matrix is singular. This is probably an issue
220  // either with the choice of rows the rows we extracted, or with the
221  // input matrix itself.
222  TEUCHOS_TEST_FOR_EXCEPTION(
223  INFO > 0, std::runtime_error, "Ifpack2::DenseContainer::factor: "
224  "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
225  "computed U factor is exactly singular. U(" << INFO << "," << INFO << ") "
226  "(one-based index i) is exactly zero. This probably means that the input "
227  "matrix has a singular diagonal block.");
228 }
229 
230 template<class MatrixType, class LocalScalarType>
231 void
233 factor ()
234 {}
235 
236 template<class MatrixType, class LocalScalarType>
237 void
239 applyImpl (const local_mv_type& X,
240  local_mv_type& Y,
241  Teuchos::ETransp mode,
242  LocalScalarType alpha,
243  LocalScalarType beta) const
244 {
245  using Teuchos::ArrayRCP;
246  using Teuchos::RCP;
247  using Teuchos::rcp;
248  using Teuchos::rcpFromRef;
249 
250  TEUCHOS_TEST_FOR_EXCEPTION(
251  X.getLocalLength () != Y.getLocalLength (),
252  std::logic_error, "Ifpack2::DenseContainer::applyImpl: X and Y have "
253  "incompatible dimensions (" << X.getLocalLength () << " resp. "
254  << Y.getLocalLength () << "). Please report this bug to "
255  "the Ifpack2 developers.");
256  TEUCHOS_TEST_FOR_EXCEPTION(
257  localMap_->getNodeNumElements () != X.getLocalLength (),
258  std::logic_error, "Ifpack2::DenseContainer::applyImpl: The inverse "
259  "operator and X have incompatible dimensions (" <<
260  localMap_->getNodeNumElements () << " resp. "
261  << X.getLocalLength () << "). Please report this bug to "
262  "the Ifpack2 developers.");
263  TEUCHOS_TEST_FOR_EXCEPTION(
264  localMap_->getNodeNumElements () != Y.getLocalLength (),
265  std::logic_error, "Ifpack2::DenseContainer::applyImpl: The inverse "
266  "operator and Y have incompatible dimensions (" <<
267  localMap_->getNodeNumElements () << " resp. "
268  << Y.getLocalLength () << "). Please report this bug to "
269  "the Ifpack2 developers.");
270  TEUCHOS_TEST_FOR_EXCEPTION(
271  X.getLocalLength () != static_cast<size_t> (diagBlock_.numRows ()),
272  std::logic_error, "Ifpack2::DenseContainer::applyImpl: The input "
273  "multivector X has incompatible dimensions from those of the "
274  "inverse operator (" << X.getLocalLength () << " vs. "
275  << (mode == Teuchos::NO_TRANS ? diagBlock_.numCols () : diagBlock_.numRows ())
276  << "). Please report this bug to the Ifpack2 developers.");
277  TEUCHOS_TEST_FOR_EXCEPTION(
278  X.getLocalLength () != static_cast<size_t> (diagBlock_.numRows ()),
279  std::logic_error, "Ifpack2::DenseContainer::applyImpl: The output "
280  "multivector Y has incompatible dimensions from those of the "
281  "inverse operator (" << Y.getLocalLength () << " vs. "
282  << (mode == Teuchos::NO_TRANS ? diagBlock_.numRows () : diagBlock_.numCols ())
283  << "). Please report this bug to the Ifpack2 developers.");
284 
285  typedef Teuchos::ScalarTraits<local_scalar_type> STS;
286  const int numVecs = static_cast<int> (X.getNumVectors ());
287  if (alpha == STS::zero ()) { // don't need to solve the linear system
288  if (beta == STS::zero ()) {
289  // Use BLAS AXPY semantics for beta == 0: overwrite, clobbering
290  // any Inf or NaN values in Y (rather than multiplying them by
291  // zero, resulting in NaN values).
292  Y.putScalar (STS::zero ());
293  }
294  else { // beta != 0
295  Y.scale (beta);
296  }
297  }
298  else { // alpha != 0; must solve the linear system
299  Teuchos::LAPACK<int, local_scalar_type> lapack;
300  // If beta is nonzero or Y is not constant stride, we have to use
301  // a temporary output multivector. It gets a (deep) copy of X,
302  // since GETRS overwrites its (multi)vector input with its output.
303  RCP<local_mv_type> Y_tmp;
304  if (beta == STS::zero () ){
305  Tpetra::deep_copy (Y, X);
306  Y_tmp = rcpFromRef (Y);
307  }
308  else {
309  Y_tmp = rcp (new local_mv_type (X, Teuchos::Copy));
310  }
311  const int Y_stride = static_cast<int> (Y_tmp->getStride ());
312  ArrayRCP<local_scalar_type> Y_view = Y_tmp->get1dViewNonConst ();
313  local_scalar_type* const Y_ptr = Y_view.getRawPtr ();
314 
315  int INFO = 0;
316  const char trans =
317  (mode == Teuchos::CONJ_TRANS ? 'C' : (mode == Teuchos::TRANS ? 'T' : 'N'));
318  lapack.GETRS (trans, diagBlock_.numRows (), numVecs,
319  diagBlock_.values (), diagBlock_.stride (),
320  ipiv_.getRawPtr (), Y_ptr, Y_stride, &INFO);
321  TEUCHOS_TEST_FOR_EXCEPTION(
322  INFO != 0, std::runtime_error, "Ifpack2::DenseContainer::applyImpl: "
323  "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
324  "failed with INFO = " << INFO << " != 0.");
325 
326  if (beta != STS::zero ()) {
327  Y.update (alpha, *Y_tmp, beta);
328  }
329  else if (! Y.isConstantStride ()) {
330  Tpetra::deep_copy (Y, *Y_tmp);
331  }
332  }
333 }
334 
335 template<class MatrixType, class LocalScalarType>
336 void
338 applyImpl (const local_mv_type& /* X */,
339  local_mv_type& /* Y */,
340  Teuchos::ETransp /* mode */,
341  LocalScalarType /* alpha */,
342  LocalScalarType /* beta */) const
343 {}
344 
345 template<class MatrixType, class LocalScalarType>
346 void
348 apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
349  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
350  Teuchos::ETransp mode,
351  scalar_type alpha,
352  scalar_type beta) const
353 {
354  using Teuchos::ArrayView;
355  using Teuchos::as;
356  using Teuchos::RCP;
357  using Teuchos::rcp;
358 
359  // The local operator might have a different Scalar type than
360  // MatrixType. This means that we might have to convert X and Y to
361  // the Tpetra::MultiVector specialization that the local operator
362  // wants. This class' X_ and Y_ internal fields are of the right
363  // type for the local operator, so we can use those as targets.
364 
365  // Tpetra::MultiVector specialization corresponding to the global operator.
366  typedef Tpetra::MultiVector<scalar_type,
368 
369  const char prefix[] = "Ifpack2::DenseContainer::weightedApply: ";
370  TEUCHOS_TEST_FOR_EXCEPTION(
371  ! IsComputed_, std::runtime_error, prefix << "You must have called the "
372  "compute() method before you may call this method. You may call "
373  "apply() as many times as you want after calling compute() once, "
374  "but you must have called compute() at least once first.");
375  const size_t numVecs = X.getNumVectors ();
376  TEUCHOS_TEST_FOR_EXCEPTION(
377  numVecs != Y.getNumVectors (), std::runtime_error,
378  prefix << "X and Y have different numbers of vectors (columns). X has "
379  << X.getNumVectors () << ", but Y has " << X.getNumVectors () << ".");
380 
381  if (numVecs == 0) {
382  return; // done! nothing to do
383  }
384 
385  // The local operator works on a permuted subset of the local parts
386  // of X and Y. The subset and permutation are defined by the index
387  // array returned by getLocalRows(). If the permutation is trivial
388  // and the subset is exactly equal to the local indices, then we
389  // could use the local parts of X and Y exactly, without needing to
390  // permute. Otherwise, we have to use temporary storage to permute
391  // X and Y. For now, we always use temporary storage.
392  //
393  // Create temporary permuted versions of the input and output.
394  // (Re)allocate X_ and/or Y_ only if necessary. We'll use them to
395  // store the permuted versions of X resp. Y. Note that X_local has
396  // the domain Map of the operator, which may be a permuted subset of
397  // the local Map corresponding to X.getMap(). Similarly, Y_local
398  // has the range Map of the operator, which may be a permuted subset
399  // of the local Map corresponding to Y.getMap(). numRows_ here
400  // gives the number of rows in the row Map of the local Inverse_
401  // operator.
402  //
403  // FIXME (mfh 20 Aug 2013) There might be an implicit assumption
404  // here that the row Map and the range Map of that operator are
405  // the same.
406  //
407  // FIXME (mfh 20 Aug 2013) This "local permutation" functionality
408  // really belongs in Tpetra.
409 
410  if (X_.is_null ()) {
411  X_ = rcp (new local_mv_type (localMap_, numVecs));
412  }
413  RCP<local_mv_type> X_local = X_;
414  TEUCHOS_TEST_FOR_EXCEPTION(
415  X_local->getLocalLength () != numRows_, std::logic_error,
416  "Ifpack2::DenseContainer::apply: "
417  "X_local has length " << X_local->getLocalLength () << ", which does "
418  "not match numRows_ = " << numRows_ << ". Please report this bug to "
419  "the Ifpack2 developers.");
420  ArrayView<const local_ordinal_type> localRows = this->getLocalRows ();
421 
423  mvgs.gather (*X_local, X, localRows);
424 
425  // We must gather the contents of the output multivector Y even on
426  // input to applyImpl(), since the inverse operator might use it as
427  // an initial guess for a linear solve. We have no way of knowing
428  // whether it does or does not.
429 
430  if (Y_.is_null ()) {
431  Y_ = rcp (new local_mv_type (localMap_, numVecs));
432  }
433  RCP<local_mv_type> Y_local = Y_;
434  TEUCHOS_TEST_FOR_EXCEPTION(
435  Y_local->getLocalLength () != numRows_, std::logic_error,
436  "Ifpack2::DenseContainer::apply: "
437  "Y_local has length " << X_local->getLocalLength () << ", which does "
438  "not match numRows_ = " << numRows_ << ". Please report this bug to "
439  "the Ifpack2 developers.");
440  mvgs.gather (*Y_local, Y, localRows);
441 
442  // Apply the local operator:
443  // Y_local := beta*Y_local + alpha*M^{-1}*X_local
444  this->applyImpl (*X_local, *Y_local, mode, as<local_scalar_type> (alpha),
445  as<local_scalar_type> (beta));
446 
447  // Scatter the permuted subset output vector Y_local back into the
448  // original output multivector Y.
449  mvgs.scatter (Y, *Y_local, localRows);
450 }
451 
452 template<class MatrixType, class LocalScalarType>
453 void
455 apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& /* X */,
456  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& /* Y */,
457  Teuchos::ETransp /* mode */,
458  scalar_type /* alpha */,
459  scalar_type /* beta */) const
460 {}
461 
462 template<class MatrixType, class LocalScalarType>
464 weightedApply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
465  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
466  const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& D,
467  Teuchos::ETransp mode,
468  scalar_type alpha,
469  scalar_type beta) const
470 {
471  using Teuchos::ArrayRCP;
472  using Teuchos::ArrayView;
473  using Teuchos::Range1D;
474  using Teuchos::RCP;
475  using Teuchos::rcp;
476  using Teuchos::rcp_const_cast;
477  using std::cerr;
478  using std::endl;
479  typedef Teuchos::ScalarTraits<scalar_type> STS;
480 
481  // The local operator template parameter might have a different
482  // Scalar type than MatrixType. This means that we might have to
483  // convert X and Y to the Tpetra::MultiVector specialization that
484  // the local operator wants. This class' X_ and Y_ internal fields
485  // are of the right type for the local operator, so we can use those
486  // as targets.
487 
488  // Tpetra::MultiVector specialization corresponding to the global operator.
489  typedef Tpetra::MultiVector<scalar_type,
491  // Tpetra::Vector specialization corresponding to the local
492  // operator. We will use this for the subset permutation of the
493  // diagonal scaling D.
494  typedef Tpetra::Vector<local_scalar_type, local_ordinal_type,
495  global_ordinal_type, node_type> local_vec_type;
496 
497  const char prefix[] = "Ifpack2::DenseContainer::weightedApply: ";
498  TEUCHOS_TEST_FOR_EXCEPTION(
499  ! IsComputed_, std::runtime_error, prefix << "You must have called the "
500  "compute() method before you may call this method. You may call "
501  "weightedApply() as many times as you want after calling compute() once, "
502  "but you must have called compute() at least once first.");
503  const size_t numVecs = X.getNumVectors ();
504  TEUCHOS_TEST_FOR_EXCEPTION(
505  numVecs != Y.getNumVectors (), std::runtime_error,
506  prefix << "X and Y have different numbers of vectors (columns). X has "
507  << X.getNumVectors () << ", but Y has " << X.getNumVectors () << ".");
508 
509  if (numVecs == 0) {
510  return; // done! nothing to do
511  }
512 
513  // The local operator works on a permuted subset of the local parts
514  // of X and Y. The subset and permutation are defined by the index
515  // array returned by getLocalRows(). If the permutation is trivial
516  // and the subset is exactly equal to the local indices, then we
517  // could use the local parts of X and Y exactly, without needing to
518  // permute. Otherwise, we have to use temporary storage to permute
519  // X and Y. For now, we always use temporary storage.
520  //
521  // Create temporary permuted versions of the input and output.
522  // (Re)allocate X_ and/or Y_ only if necessary. We'll use them to
523  // store the permuted versions of X resp. Y. Note that X_local has
524  // the domain Map of the operator, which may be a permuted subset of
525  // the local Map corresponding to X.getMap(). Similarly, Y_local
526  // has the range Map of the operator, which may be a permuted subset
527  // of the local Map corresponding to Y.getMap(). numRows_ here
528  // gives the number of rows in the row Map of the local operator.
529  //
530  // FIXME (mfh 20 Aug 2013) There might be an implicit assumption
531  // here that the row Map and the range Map of that operator are
532  // the same.
533  //
534  // FIXME (mfh 20 Aug 2013) This "local permutation" functionality
535  // really belongs in Tpetra.
536 
537  if (X_.is_null ()) {
538  X_ = rcp (new local_mv_type (localMap_, numVecs));
539  }
540  RCP<local_mv_type> X_local = X_;
541  TEUCHOS_TEST_FOR_EXCEPTION(
542  X_local->getLocalLength () != numRows_, std::logic_error,
543  "Ifpack2::DenseContainer::weightedApply: "
544  "X_local has length " << X_local->getLocalLength () << ", which does "
545  "not match numRows_ = " << numRows_ << ". Please report this bug to "
546  "the Ifpack2 developers.");
547  ArrayView<const local_ordinal_type> localRows = this->getLocalRows ();
548 
550  mvgs.gather (*X_local, X, localRows);
551 
552  // We must gather the output multivector Y even on input to
553  // applyImpl(), since the local operator might use it as an initial
554  // guess for a linear solve. We have no way of knowing whether it
555  // does or does not.
556 
557  if (Y_.is_null ()) {
558  Y_ = rcp (new local_mv_type (localMap_, numVecs));
559  }
560  RCP<local_mv_type> Y_local = Y_;
561  TEUCHOS_TEST_FOR_EXCEPTION(
562  Y_local->getLocalLength () != numRows_, std::logic_error,
563  "Ifpack2::DenseContainer::weightedApply: "
564  "Y_local has length " << X_local->getLocalLength () << ", which does "
565  "not match numRows_ = " << numRows_ << ". Please report this bug to "
566  "the Ifpack2 developers.");
567  mvgs.gather (*Y_local, Y, localRows);
568 
569  // Apply the diagonal scaling D to the input X. It's our choice
570  // whether the result has the original input Map of X, or the
571  // permuted subset Map of X_local. If the latter, we also need to
572  // gather D into the permuted subset Map. We choose the latter, to
573  // save memory and computation. Thus, we do the following:
574  //
575  // 1. Gather D into a temporary vector D_local.
576  // 2. Create a temporary X_scaled to hold diag(D_local) * X_local.
577  // 3. Compute X_scaled := diag(D_loca) * X_local.
578 
579  local_vec_type D_local (localMap_);
580  TEUCHOS_TEST_FOR_EXCEPTION(
581  D_local.getLocalLength () != numRows_, std::logic_error,
582  "Ifpack2::DenseContainer::weightedApply: "
583  "D_local has length " << D_local.getLocalLength () << ", which does "
584  "not match numRows_ = " << numRows_ << ". Please report this bug to "
585  "the Ifpack2 developers.");
586  mvgs.gather (D_local, D, localRows);
587  local_mv_type X_scaled (localMap_, numVecs);
588  X_scaled.elementWiseMultiply (STS::one (), D_local, *X_local, STS::zero ());
589 
590  // Y_temp will hold the result of M^{-1}*X_scaled. If beta == 0, we
591  // can write the result of Inverse_->apply() directly to Y_local, so
592  // Y_temp may alias Y_local. Otherwise, if beta != 0, we need
593  // temporary storage for M^{-1}*X_scaled, so Y_temp must be
594  // different than Y_local.
595  RCP<local_mv_type> Y_temp;
596  if (beta == STS::zero ()) {
597  Y_temp = Y_local;
598  } else {
599  Y_temp = rcp (new local_mv_type (localMap_, numVecs));
600  }
601 
602  // Apply the local operator: Y_temp := M^{-1} * X_scaled
603  this->applyImpl (X_scaled, *Y_temp, mode, STS::one (), STS::zero ());
604  // Y_local := beta * Y_local + alpha * diag(D_local) * Y_temp.
605  //
606  // Note that we still use the permuted subset scaling D_local here,
607  // because Y_temp has the same permuted subset Map. That's good, in
608  // fact, because it's a subset: less data to read and multiply.
609  Y_local->elementWiseMultiply (alpha, D_local, *Y_temp, beta);
610 
611  // Copy the permuted subset output vector Y_local into the original
612  // output multivector Y.
613  mvgs.scatter (Y, *Y_local, localRows);
614 }
615 
616 template<class MatrixType, class LocalScalarType>
617 void
619 weightedApply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& /* X */,
620  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& /* Y */,
621  const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& /* D */,
622  Teuchos::ETransp /* mode */,
623  scalar_type /* alpha */,
624  scalar_type /* beta */ ) const
625 {
626 }
627 
628 template<class MatrixType, class LocalScalarType>
629 std::ostream&
631 print (std::ostream& os) const
632 {
633  Teuchos::FancyOStream fos (Teuchos::rcpFromRef (os));
634  fos.setOutputToRootOnly (0);
635  this->describe (fos);
636  return os;
637 }
638 
639 template<class MatrixType, class LocalScalarType>
640 std::ostream&
642 print (std::ostream& os) const
643 {
644  return os;
645 }
646 
647 template<class MatrixType, class LocalScalarType>
648 std::string
650 description () const
651 {
652  std::ostringstream oss;
653  oss << Teuchos::Describable::description();
654  if (isInitialized()) {
655  if (isComputed()) {
656  oss << "{status = initialized, computed";
657  }
658  else {
659  oss << "{status = initialized, not computed";
660  }
661  }
662  else {
663  oss << "{status = not initialized, not computed";
664  }
665 
666  oss << "}";
667  return oss.str();
668 }
669 
670 template<class MatrixType, class LocalScalarType>
671 std::string
673 description () const
674 {
675  return "";
676 }
677 
678 template<class MatrixType, class LocalScalarType>
679 void
681 describe (Teuchos::FancyOStream& os,
682  const Teuchos::EVerbosityLevel verbLevel) const
683 {
684  using std::endl;
685  if(verbLevel==Teuchos::VERB_NONE) return;
686  os << "================================================================================" << endl;
687  os << "Ifpack2::DenseContainer" << endl;
688  os << "Number of rows = " << numRows_ << endl;
689  os << "isInitialized() = " << IsInitialized_ << endl;
690  os << "isComputed() = " << IsComputed_ << endl;
691  os << "================================================================================" << endl;
692  os << endl;
693 }
694 
695 template<class MatrixType, class LocalScalarType>
696 void
698 describe (Teuchos::FancyOStream& os,
699  const Teuchos::EVerbosityLevel verbLevel) const
700 {}
701 
702 template<class MatrixType, class LocalScalarType>
703 void
705 extract (const Teuchos::RCP<const row_matrix_type>& globalMatrix)
706 {
707  using Teuchos::Array;
708  using Teuchos::ArrayView;
709  using Teuchos::toString;
710  typedef local_ordinal_type LO;
711  typedef global_ordinal_type GO;
712  typedef Tpetra::Map<LO, GO, node_type> map_type;
713  const size_t inputMatrixNumRows = globalMatrix->getNodeNumRows ();
714  // We only use the rank of the calling process and the number of MPI
715  // processes for generating error messages. Extraction itself is
716  // entirely local to each participating MPI process.
717  const int myRank = globalMatrix->getRowMap ()->getComm ()->getRank ();
718  const int numProcs = globalMatrix->getRowMap ()->getComm ()->getSize ();
719 
720  // Sanity check that the local row indices to extract fall within
721  // the valid range of local row indices for the input matrix.
722  ArrayView<const LO> localRows = this->getLocalRows ();
723  for (size_t j = 0; j < numRows_; ++j) {
724  TEUCHOS_TEST_FOR_EXCEPTION(
725  localRows[j] < 0 ||
726  static_cast<size_t> (localRows[j]) >= inputMatrixNumRows,
727  std::runtime_error, "Ifpack2::DenseContainer::extract: On process " <<
728  myRank << " of " << numProcs << ", localRows[j=" << j << "] = " <<
729  localRows[j] << ", which is out of the valid range of local row indices "
730  "indices [0, " << (inputMatrixNumRows - 1) << "] for the input matrix.");
731  }
732 
733  // Convert the local row indices we want into local column indices.
734  // For every local row ii_local = localRows[i] we take, we also want
735  // to take the corresponding column. To find the corresponding
736  // column, we use the row Map to convert the local row index
737  // ii_local into a global index ii_global, and then use the column
738  // Map to convert ii_global into a local column index jj_local. If
739  // the input matrix doesn't have a column Map, we need to be using
740  // global indices anyway...
741 
742  // We use the domain Map to exclude off-process global entries.
743  const map_type& globalRowMap = * (globalMatrix->getRowMap ());
744  const map_type& globalColMap = * (globalMatrix->getColMap ());
745  const map_type& globalDomMap = * (globalMatrix->getDomainMap ());
746 
747  bool rowIndsValid = true;
748  bool colIndsValid = true;
749  Array<LO> localCols (numRows_);
750  // For error messages, collect the sets of invalid row indices and
751  // invalid column indices. They are otherwise not useful.
752  Array<LO> invalidLocalRowInds;
753  Array<GO> invalidGlobalColInds;
754  for (size_t i = 0; i < numRows_; ++i) {
755  // ii_local is the (local) row index we want to look up.
756  const LO ii_local = localRows[i];
757  // Find the global index jj_global corresponding to ii_local.
758  // Global indices are the same (rather, are required to be the
759  // same) in all three Maps, which is why we use jj (suggesting a
760  // column index, which is how we will use it below).
761  const GO jj_global = globalRowMap.getGlobalElement (ii_local);
762  if (jj_global == Teuchos::OrdinalTraits<GO>::invalid ()) {
763  // If ii_local is not a local index in the row Map on the
764  // calling process, that means localRows is incorrect. We've
765  // already checked for this in the constructor, but we might as
766  // well check again here, since it's cheap to do so (just an
767  // integer comparison, since we need jj_global anyway).
768  rowIndsValid = false;
769  invalidLocalRowInds.push_back (ii_local);
770  break;
771  }
772  // Exclude "off-process" entries: that is, those in the column Map
773  // on this process that are not in the domain Map on this process.
774  if (globalDomMap.isNodeGlobalElement (jj_global)) {
775  // jj_global is not an off-process entry. Look up its local
776  // index in the column Map; we want to extract this column index
777  // from the input matrix. If jj_global is _not_ in the column
778  // Map on the calling process, that could mean that the column
779  // in question is empty on this process. That would be bad for
780  // solving linear systems with the extract submatrix. We could
781  // solve the resulting singular linear systems in a minimum-norm
782  // least-squares sense, but for now we simply raise an exception.
783  const LO jj_local = globalColMap.getLocalElement (jj_global);
784  if (jj_local == Teuchos::OrdinalTraits<local_ordinal_type>::invalid ()) {
785  colIndsValid = false;
786  invalidGlobalColInds.push_back (jj_global);
787  break;
788  }
789  localCols[i] = jj_local;
790  }
791  }
792  TEUCHOS_TEST_FOR_EXCEPTION(
793  ! rowIndsValid, std::logic_error, "Ifpack2::DenseContainer::extract: "
794  "On process " << myRank << ", at least one row index in the set of local "
795  "row indices given to the constructor is not a valid local row index in "
796  "the input matrix's row Map on this process. This should be impossible "
797  "because the constructor checks for this case. Here is the complete set "
798  "of invalid local row indices: " << toString (invalidLocalRowInds) << ". "
799  "Please report this bug to the Ifpack2 developers.");
800  TEUCHOS_TEST_FOR_EXCEPTION(
801  ! colIndsValid, std::runtime_error, "Ifpack2::DenseContainer::extract: "
802  "On process " << myRank << ", "
803  "At least one row index in the set of row indices given to the constructor "
804  "does not have a corresponding column index in the input matrix's column "
805  "Map. This probably means that the column(s) in question is/are empty on "
806  "this process, which would make the submatrix to extract structurally "
807  "singular. Here is the compete set of invalid global column indices: "
808  << toString (invalidGlobalColInds) << ".");
809 
810  diagBlock_.putScalar (Teuchos::ScalarTraits<local_scalar_type>::zero ());
811 
812  const size_t maxNumEntriesInRow = globalMatrix->getNodeMaxNumRowEntries ();
813  Array<scalar_type> val (maxNumEntriesInRow);
814  Array<local_ordinal_type> ind (maxNumEntriesInRow);
815 
816  const local_ordinal_type INVALID =
817  Teuchos::OrdinalTraits<local_ordinal_type>::invalid ();
818  for (size_t i = 0; i < numRows_; ++i) {
819  const local_ordinal_type localRow = localRows[i];
820  size_t numEntries;
821  globalMatrix->getLocalRowCopy (localRow, ind (), val (), numEntries);
822 
823  for (size_t k = 0; k < numEntries; ++k) {
824  const local_ordinal_type localCol = ind[k];
825 
826  // Skip off-process elements
827  //
828  // FIXME (mfh 24 Aug 2013) This assumes the following:
829  //
830  // 1. The column and row Maps begin with the same set of
831  // on-process entries, in the same order. That is,
832  // on-process row and column indices are the same.
833  // 2. All off-process indices in the column Map of the input
834  // matrix occur after that initial set.
835  if (localCol >= 0 && static_cast<size_t> (localCol) < inputMatrixNumRows) {
836  // for local column IDs, look for each ID in the list
837  // of columns hosted by this object
838  local_ordinal_type jj = INVALID;
839  for (size_t kk = 0; kk < numRows_; ++kk) {
840  if (localRows[kk] == localCol) {
841  jj = kk;
842  }
843  }
844  if (jj != INVALID) {
845  diagBlock_ (i, jj) += val[k]; // ???
846  }
847  }
848  }
849  }
850 }
851 
852 template<class MatrixType, class LocalScalarType>
853 void
855 extract (const Teuchos::RCP<const row_matrix_type>& /* globalMatrix */)
856 {}
857 
858 } // namespace Ifpack2
859 
860 // There's no need to instantiate for CrsMatrix too. All Ifpack2
861 // preconditioners can and should do dynamic casts if they need a type
862 // more specific than RowMatrix.
863 
864 #define IFPACK2_DENSECONTAINER_INSTANT(S,LO,GO,N) \
865  template class Ifpack2::DenseContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
866 
867 #endif // IFPACK2_SPARSECONTAINER_HPP
DenseContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::ArrayView< const local_ordinal_type > &localRows)
Constructor.
Definition: Ifpack2_DenseContainer_def.hpp:61
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input (global) matrix.
Definition: Ifpack2_DenseContainer_decl.hpp:332
Store and solve a local dense linear problem.
Definition: Ifpack2_DenseContainer_decl.hpp:108
LocalScalarType local_scalar_type
The second template parameter of this class.
Definition: Ifpack2_DenseContainer_decl.hpp:325
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input (global) matrix.
Definition: Ifpack2_DenseContainer_decl.hpp:134
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input (global) matrix.
Definition: Ifpack2_DenseContainer_decl.hpp:330
LocalScalarType local_scalar_type
The second template parameter of this class.
Definition: Ifpack2_DenseContainer_decl.hpp:127
MatrixType::node_type node_type
The Node type of the input (global) matrix.
Definition: Ifpack2_DenseContainer_decl.hpp:136
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input (global) matrix.
Definition: Ifpack2_DenseContainer_decl.hpp:132
Interface for creating and solving a local linear problem.
Definition: Ifpack2_Container.hpp:103
BELOS_DEPRECATED const char * toString(const StatusType status)
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
Implementation detail of Ifpack2::Container subclasses.
Definition: Ifpack2_Details_MultiVectorLocalGatherScatter.hpp:85