Ifpack2 Templated Preconditioning Package  Version 1.0
Ifpack2_ILUT_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_ILUT_DEF_HPP
44 #define IFPACK2_ILUT_DEF_HPP
45 
46 // disable clang warnings
47 #if defined (__clang__) && !defined (__INTEL_COMPILER)
48 #pragma clang system_header
49 #endif
50 
51 #include <Ifpack2_Heap.hpp>
52 #include <Ifpack2_LocalFilter.hpp>
53 #include <Ifpack2_Parameters.hpp>
54 #include <Tpetra_CrsMatrix.hpp>
55 #include <Teuchos_Time.hpp>
56 #include <Teuchos_TypeNameTraits.hpp>
57 
58 
59 namespace Ifpack2 {
60 
61  namespace {
62 
87  template<class ScalarType>
88  typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
89  ilutDefaultDropTolerance () {
90  using Teuchos::as;
91  typedef Teuchos::ScalarTraits<ScalarType> STS;
92  typedef typename STS::magnitudeType magnitude_type;
93  typedef Teuchos::ScalarTraits<magnitude_type> STM;
94 
95  // 1/2. Hopefully this can be represented in magnitude_type.
96  const magnitude_type oneHalf = STM::one() / (STM::one() + STM::one());
97 
98  // The min ensures that in case magnitude_type has very low
99  // precision, we'll at least get some value strictly less than
100  // one.
101  return std::min (as<magnitude_type> (1000) * STS::magnitude (STS::eps ()), oneHalf);
102  }
103 
104  // Full specialization for ScalarType = double.
105  // This specialization preserves ILUT's previous default behavior.
106  template<>
107  Teuchos::ScalarTraits<double>::magnitudeType
108  ilutDefaultDropTolerance<double> () {
109  return 1e-12;
110  }
111 
112  } // namespace (anonymous)
113 
114 
115 template <class MatrixType>
116 ILUT<MatrixType>::ILUT (const Teuchos::RCP<const row_matrix_type>& A) :
117  A_ (A),
118  Athresh_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
119  Rthresh_ (Teuchos::ScalarTraits<magnitude_type>::one ()),
120  RelaxValue_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
121  LevelOfFill_ (1),
122  DropTolerance_ (ilutDefaultDropTolerance<scalar_type> ()),
123  InitializeTime_ (0.0),
124  ComputeTime_ (0.0),
125  ApplyTime_ (0.0),
126  NumInitialize_ (0),
127  NumCompute_ (0),
128  NumApply_ (0),
129  IsInitialized_ (false),
130  IsComputed_ (false)
131 {}
132 
133 template <class MatrixType>
135 {}
136 
137 template <class MatrixType>
138 void ILUT<MatrixType>::setParameters (const Teuchos::ParameterList& params)
139 {
140  using Teuchos::as;
141  using Teuchos::Exceptions::InvalidParameterName;
142  using Teuchos::Exceptions::InvalidParameterType;
143 
144  // Default values of the various parameters.
145  int fillLevel = 1;
146  magnitude_type absThresh = STM::zero ();
147  magnitude_type relThresh = STM::one ();
148  magnitude_type relaxValue = STM::zero ();
149  magnitude_type dropTol = ilutDefaultDropTolerance<scalar_type> ();
150 
151  bool gotFillLevel = false;
152  try {
153  // Try getting the fill level as an int.
154  fillLevel = params.get<int> ("fact: ilut level-of-fill");
155  gotFillLevel = true;
156  }
157  catch (InvalidParameterName&) {
158  // We didn't really get it, but this tells us to stop looking.
159  gotFillLevel = true;
160  }
161  catch (InvalidParameterType&) {
162  // The name is right, but the type is wrong; try different types.
163  // We don't have to check InvalidParameterName again, since we
164  // checked it above, and it has nothing to do with the type.
165  }
166 
167  if (! gotFillLevel) {
168  // Try magnitude_type, for backwards compatibility.
169  try {
170  fillLevel = as<int> (params.get<magnitude_type> ("fact: ilut level-of-fill"));
171  }
172  catch (InvalidParameterType&) {}
173  }
174  if (! gotFillLevel) {
175  // Try double, for backwards compatibility.
176  try {
177  fillLevel = as<int> (params.get<double> ("fact: ilut level-of-fill"));
178  }
179  catch (InvalidParameterType&) {}
180  }
181  // If none of the above attempts succeed, accept the default value.
182 
183  TEUCHOS_TEST_FOR_EXCEPTION(
184  fillLevel <= 0, std::runtime_error,
185  "Ifpack2::ILUT: The \"fact: ilut level-of-fill\" parameter must be "
186  "strictly greater than zero, but you specified a value of " << fillLevel
187  << ". Remember that for ILUT, the fill level p means something different "
188  "than it does for ILU(k). ILU(0) produces factors with the same sparsity "
189  "structure as the input matrix A; ILUT with p = 0 always produces a "
190  "diagonal matrix, and is thus probably not what you want.");
191 
192  try {
193  absThresh = params.get<magnitude_type> ("fact: absolute threshold");
194  }
195  catch (InvalidParameterType&) {
196  // Try double, for backwards compatibility.
197  // The cast from double to magnitude_type must succeed.
198  absThresh = as<magnitude_type> (params.get<double> ("fact: absolute threshold"));
199  }
200  catch (InvalidParameterName&) {
201  // Accept the default value.
202  }
203 
204  try {
205  relThresh = params.get<magnitude_type> ("fact: relative threshold");
206  }
207  catch (InvalidParameterType&) {
208  // Try double, for backwards compatibility.
209  // The cast from double to magnitude_type must succeed.
210  relThresh = as<magnitude_type> (params.get<double> ("fact: relative threshold"));
211  }
212  catch (InvalidParameterName&) {
213  // Accept the default value.
214  }
215 
216  try {
217  relaxValue = params.get<magnitude_type> ("fact: relax value");
218  }
219  catch (InvalidParameterType&) {
220  // Try double, for backwards compatibility.
221  // The cast from double to magnitude_type must succeed.
222  relaxValue = as<magnitude_type> (params.get<double> ("fact: relax value"));
223  }
224  catch (InvalidParameterName&) {
225  // Accept the default value.
226  }
227 
228  try {
229  dropTol = params.get<magnitude_type> ("fact: drop tolerance");
230  }
231  catch (InvalidParameterType&) {
232  // Try double, for backwards compatibility.
233  // The cast from double to magnitude_type must succeed.
234  dropTol = as<magnitude_type> (params.get<double> ("fact: drop tolerance"));
235  }
236  catch (InvalidParameterName&) {
237  // Accept the default value.
238  }
239 
240  // "Commit" the values only after validating all of them. This
241  // ensures that there are no side effects if this routine throws an
242  // exception.
243 
244  // mfh 28 Nov 2012: The previous code would not assign Athresh_,
245  // Rthresh_, RelaxValue_, or DropTolerance_ if the read-in value was
246  // -1. I don't know if keeping this behavior is correct, but I'll
247  // keep it just so as not to change previous behavior.
248 
249  LevelOfFill_ = fillLevel;
250  if (absThresh != -STM::one ()) {
251  Athresh_ = absThresh;
252  }
253  if (relThresh != -STM::one ()) {
254  Rthresh_ = relThresh;
255  }
256  if (relaxValue != -STM::one ()) {
257  RelaxValue_ = relaxValue;
258  }
259  if (dropTol != -STM::one ()) {
260  DropTolerance_ = dropTol;
261  }
262 }
263 
264 
265 template <class MatrixType>
266 Teuchos::RCP<const Teuchos::Comm<int> >
268  TEUCHOS_TEST_FOR_EXCEPTION(
269  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getComm: "
270  "The matrix is null. Please call setMatrix() with a nonnull input "
271  "before calling this method.");
272  return A_->getComm ();
273 }
274 
275 
276 template <class MatrixType>
277 Teuchos::RCP<const typename ILUT<MatrixType>::row_matrix_type>
279  return A_;
280 }
281 
282 
283 template <class MatrixType>
284 Teuchos::RCP<const typename ILUT<MatrixType>::map_type>
286 {
287  TEUCHOS_TEST_FOR_EXCEPTION(
288  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getDomainMap: "
289  "The matrix is null. Please call setMatrix() with a nonnull input "
290  "before calling this method.");
291  return A_->getDomainMap ();
292 }
293 
294 
295 template <class MatrixType>
296 Teuchos::RCP<const typename ILUT<MatrixType>::map_type>
298 {
299  TEUCHOS_TEST_FOR_EXCEPTION(
300  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getRangeMap: "
301  "The matrix is null. Please call setMatrix() with a nonnull input "
302  "before calling this method.");
303  return A_->getRangeMap ();
304 }
305 
306 
307 template <class MatrixType>
309  return true;
310 }
311 
312 
313 template <class MatrixType>
315  return NumInitialize_;
316 }
317 
318 
319 template <class MatrixType>
321  return NumCompute_;
322 }
323 
324 
325 template <class MatrixType>
327  return NumApply_;
328 }
329 
330 
331 template <class MatrixType>
333  return InitializeTime_;
334 }
335 
336 
337 template<class MatrixType>
339  return ComputeTime_;
340 }
341 
342 
343 template<class MatrixType>
345  return ApplyTime_;
346 }
347 
348 
349 template<class MatrixType>
351  return L_->getGlobalNumEntries () + U_->getGlobalNumEntries ();
352 }
353 
354 
355 template<class MatrixType>
357  return L_->getNodeNumEntries () + U_->getNodeNumEntries ();
358 }
359 
360 
361 template<class MatrixType>
362 void ILUT<MatrixType>::setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
363 {
364  if (A.getRawPtr () != A_.getRawPtr ()) {
365  // Check in serial or one-process mode if the matrix is square.
366  TEUCHOS_TEST_FOR_EXCEPTION(
367  ! A.is_null () && A->getComm ()->getSize () == 1 &&
368  A->getNodeNumRows () != A->getNodeNumCols (),
369  std::runtime_error, "Ifpack2::ILUT::setMatrix: If A's communicator only "
370  "contains one process, then A must be square. Instead, you provided a "
371  "matrix A with " << A->getNodeNumRows () << " rows and "
372  << A->getNodeNumCols () << " columns.");
373 
374  // It's legal for A to be null; in that case, you may not call
375  // initialize() until calling setMatrix() with a nonnull input.
376  // Regardless, setting the matrix invalidates any previous
377  // factorization.
378  IsInitialized_ = false;
379  IsComputed_ = false;
380  A_local_ = Teuchos::null;
381  L_ = Teuchos::null;
382  U_ = Teuchos::null;
383  A_ = A;
384  }
385 }
386 
387 
388 template<class MatrixType>
390 {
391  Teuchos::Time timer ("ILUT::initialize");
392  {
393  Teuchos::TimeMonitor timeMon (timer);
394 
395  // Check that the matrix is nonnull.
396  TEUCHOS_TEST_FOR_EXCEPTION(
397  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::initialize: "
398  "The matrix to precondition is null. Please call setMatrix() with a "
399  "nonnull input before calling this method.");
400 
401  // Clear any previous computations.
402  IsInitialized_ = false;
403  IsComputed_ = false;
404  A_local_ = Teuchos::null;
405  L_ = Teuchos::null;
406  U_ = Teuchos::null;
407 
408  A_local_ = makeLocalFilter (A_); // Compute the local filter.
409 
410  IsInitialized_ = true;
411  ++NumInitialize_;
412  }
413  InitializeTime_ += timer.totalElapsedTime ();
414 }
415 
416 
417 template<typename ScalarType>
418 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
419 scalar_mag (const ScalarType& s)
420 {
421  return Teuchos::ScalarTraits<ScalarType>::magnitude(s);
422 }
423 
424 
425 template<class MatrixType>
427 {
428  using Teuchos::Array;
429  using Teuchos::ArrayRCP;
430  using Teuchos::ArrayView;
431  using Teuchos::as;
432  using Teuchos::rcp;
433  using Teuchos::reduceAll;
434 
435  //--------------------------------------------------------------------------
436  // Ifpack2::ILUT is a translation of the Aztec ILUT implementation. The Aztec
437  // ILUT implementation was written by Ray Tuminaro.
438  //
439  // This isn't an exact translation of the Aztec ILUT algorithm, for the
440  // following reasons:
441  // 1. Minor differences result from the fact that Aztec factors a MSR format
442  // matrix in place, while the code below factors an input CrsMatrix which
443  // remains untouched and stores the resulting factors in separate L and U
444  // CrsMatrix objects.
445  // Also, the Aztec code begins by shifting the matrix pointers back
446  // by one, and the pointer contents back by one, and then using 1-based
447  // Fortran-style indexing in the algorithm. This Ifpack2 code uses C-style
448  // 0-based indexing throughout.
449  // 2. Aztec stores the inverse of the diagonal of U. This Ifpack2 code
450  // stores the non-inverted diagonal in U.
451  // The triangular solves (in Ifpack2::ILUT::apply()) are performed by
452  // calling the Tpetra::CrsMatrix::solve method on the L and U objects, and
453  // this requires U to contain the non-inverted diagonal.
454  //
455  // ABW.
456  //--------------------------------------------------------------------------
457 
458  // Don't count initialization in the compute() time.
459  if (! isInitialized ()) {
460  initialize ();
461  }
462 
463  Teuchos::Time timer ("ILUT::compute");
464  { // Timer scope for timing compute()
465  Teuchos::TimeMonitor timeMon (timer, true);
466  const scalar_type zero = STS::zero ();
467  const scalar_type one = STS::one ();
468 
469  const local_ordinal_type myNumRows = A_local_->getNodeNumRows ();
470  L_ = rcp (new crs_matrix_type (A_local_->getRowMap (), A_local_->getColMap (), 0));
471  U_ = rcp (new crs_matrix_type (A_local_->getRowMap (), A_local_->getColMap (), 0));
472 
473  // CGB: note, this caching approach may not be necessary anymore
474  // We will store ArrayView objects that are views of the rows of U, so that
475  // we don't have to repeatedly retrieve the view for each row. These will
476  // be populated row by row as the factorization proceeds.
477  Array<ArrayView<const local_ordinal_type> > Uindices (myNumRows);
478  Array<ArrayView<const scalar_type> > Ucoefs (myNumRows);
479 
480  // If this macro is defined, files containing the L and U factors
481  // will be written. DON'T CHECK IN THE CODE WITH THIS MACRO ENABLED!!!
482  // #define IFPACK2_WRITE_FACTORS
483 #ifdef IFPACK2_WRITE_FACTORS
484  std::ofstream ofsL("L.tif.mtx", std::ios::out);
485  std::ofstream ofsU("U.tif.mtx", std::ios::out);
486 #endif
487 
488  // The code here uses double for fill calculations, even though
489  // the fill level is actually an integer. The point is to avoid
490  // rounding and overflow for integer calculations. If int is <=
491  // 32 bits, it can never overflow double, so this cast is always
492  // OK as long as int has <= 32 bits.
493 
494  // Calculate how much fill will be allowed in addition to the
495  // space that corresponds to the input matrix entries.
496  double local_nnz = static_cast<double> (A_local_->getNodeNumEntries ());
497  double fill;
498  {
499  const double fillLevel = as<double> (getLevelOfFill ());
500  fill = ((fillLevel - 1) * local_nnz) / (2 * myNumRows);
501  }
502 
503  // std::ceil gives the smallest integer larger than the argument.
504  // this may give a slightly different result than Aztec's fill value in
505  // some cases.
506  double fill_ceil=std::ceil(fill);
507 
508  // Similarly to Aztec, we will allow the same amount of fill for each
509  // row, half in L and half in U.
510  size_type fillL = static_cast<size_type>(fill_ceil);
511  size_type fillU = static_cast<size_type>(fill_ceil);
512 
513  Array<scalar_type> InvDiagU (myNumRows, zero);
514 
515  Array<local_ordinal_type> tmp_idx;
516  Array<scalar_type> tmpv;
517 
518  enum { UNUSED, ORIG, FILL };
519  local_ordinal_type max_col = myNumRows;
520 
521  Array<int> pattern(max_col, UNUSED);
522  Array<scalar_type> cur_row(max_col, zero);
523  Array<magnitude_type> unorm(max_col);
524  magnitude_type rownorm;
525  Array<local_ordinal_type> L_cols_heap;
526  Array<local_ordinal_type> U_cols;
527  Array<local_ordinal_type> L_vals_heap;
528  Array<local_ordinal_type> U_vals_heap;
529 
530  // A comparison object which will be used to create 'heaps' of indices
531  // that are ordered according to the corresponding values in the
532  // 'cur_row' array.
533  greater_indirect<scalar_type,local_ordinal_type> vals_comp(cur_row);
534 
535  // =================== //
536  // start factorization //
537  // =================== //
538 
539  ArrayRCP<local_ordinal_type> ColIndicesARCP;
540  ArrayRCP<scalar_type> ColValuesARCP;
541  if (! A_local_->supportsRowViews ()) {
542  const size_t maxnz = A_local_->getNodeMaxNumRowEntries ();
543  ColIndicesARCP.resize (maxnz);
544  ColValuesARCP.resize (maxnz);
545  }
546 
547  for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
548  ArrayView<const local_ordinal_type> ColIndicesA;
549  ArrayView<const scalar_type> ColValuesA;
550  size_t RowNnz;
551 
552  if (A_local_->supportsRowViews ()) {
553  A_local_->getLocalRowView (row_i, ColIndicesA, ColValuesA);
554  RowNnz = ColIndicesA.size ();
555  }
556  else {
557  A_local_->getLocalRowCopy (row_i, ColIndicesARCP (), ColValuesARCP (), RowNnz);
558  ColIndicesA = ColIndicesARCP (0, RowNnz);
559  ColValuesA = ColValuesARCP (0, RowNnz);
560  }
561 
562  // Always include the diagonal in the U factor. The value should get
563  // set in the next loop below.
564  U_cols.push_back(row_i);
565  cur_row[row_i] = zero;
566  pattern[row_i] = ORIG;
567 
568  size_type L_cols_heaplen = 0;
569  rownorm = STM::zero ();
570  for (size_t i = 0; i < RowNnz; ++i) {
571  if (ColIndicesA[i] < myNumRows) {
572  if (ColIndicesA[i] < row_i) {
573  add_to_heap(ColIndicesA[i], L_cols_heap, L_cols_heaplen);
574  }
575  else if (ColIndicesA[i] > row_i) {
576  U_cols.push_back(ColIndicesA[i]);
577  }
578 
579  cur_row[ColIndicesA[i]] = ColValuesA[i];
580  pattern[ColIndicesA[i]] = ORIG;
581  rownorm += scalar_mag(ColValuesA[i]);
582  }
583  }
584 
585  // Alter the diagonal according to the absolute-threshold and
586  // relative-threshold values. If not set, those values default
587  // to zero and one respectively.
588  const magnitude_type rthresh = getRelativeThreshold();
589  const scalar_type& v = cur_row[row_i];
590  cur_row[row_i] = as<scalar_type> (getAbsoluteThreshold() * IFPACK2_SGN(v)) + rthresh*v;
591 
592  size_type orig_U_len = U_cols.size();
593  RowNnz = L_cols_heap.size() + orig_U_len;
594  rownorm = getDropTolerance() * rownorm/RowNnz;
595 
596  // The following while loop corresponds to the 'L30' goto's in Aztec.
597  size_type L_vals_heaplen = 0;
598  while (L_cols_heaplen > 0) {
599  local_ordinal_type row_k = L_cols_heap.front();
600 
601  scalar_type multiplier = cur_row[row_k] * InvDiagU[row_k];
602  cur_row[row_k] = multiplier;
603  magnitude_type mag_mult = scalar_mag(multiplier);
604  if (mag_mult*unorm[row_k] < rownorm) {
605  pattern[row_k] = UNUSED;
606  rm_heap_root(L_cols_heap, L_cols_heaplen);
607  continue;
608  }
609  if (pattern[row_k] != ORIG) {
610  if (L_vals_heaplen < fillL) {
611  add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
612  }
613  else if (L_vals_heaplen==0 ||
614  mag_mult < scalar_mag(cur_row[L_vals_heap.front()])) {
615  pattern[row_k] = UNUSED;
616  rm_heap_root(L_cols_heap, L_cols_heaplen);
617  continue;
618  }
619  else {
620  pattern[L_vals_heap.front()] = UNUSED;
621  rm_heap_root(L_vals_heap, L_vals_heaplen, vals_comp);
622  add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
623  }
624  }
625 
626  /* Reduce current row */
627 
628  ArrayView<const local_ordinal_type>& ColIndicesU = Uindices[row_k];
629  ArrayView<const scalar_type>& ColValuesU = Ucoefs[row_k];
630  size_type ColNnzU = ColIndicesU.size();
631 
632  for(size_type j=0; j<ColNnzU; ++j) {
633  if (ColIndicesU[j] > row_k) {
634  scalar_type tmp = multiplier * ColValuesU[j];
635  local_ordinal_type col_j = ColIndicesU[j];
636  if (pattern[col_j] != UNUSED) {
637  cur_row[col_j] -= tmp;
638  }
639  else if (scalar_mag(tmp) > rownorm) {
640  cur_row[col_j] = -tmp;
641  pattern[col_j] = FILL;
642  if (col_j > row_i) {
643  U_cols.push_back(col_j);
644  }
645  else {
646  add_to_heap(col_j, L_cols_heap, L_cols_heaplen);
647  }
648  }
649  }
650  }
651 
652  rm_heap_root(L_cols_heap, L_cols_heaplen);
653  }//end of while(L_cols_heaplen) loop
654 
655 
656  // Put indices and values for L into arrays and then into the L_ matrix.
657 
658  // first, the original entries from the L section of A:
659  for (size_type i = 0; i < ColIndicesA.size (); ++i) {
660  if (ColIndicesA[i] < row_i) {
661  tmp_idx.push_back(ColIndicesA[i]);
662  tmpv.push_back(cur_row[ColIndicesA[i]]);
663  pattern[ColIndicesA[i]] = UNUSED;
664  }
665  }
666 
667  // next, the L entries resulting from fill:
668  for (size_type j = 0; j < L_vals_heaplen; ++j) {
669  tmp_idx.push_back(L_vals_heap[j]);
670  tmpv.push_back(cur_row[L_vals_heap[j]]);
671  pattern[L_vals_heap[j]] = UNUSED;
672  }
673 
674  // L has a one on the diagonal, but we don't explicitly store
675  // it. If we don't store it, then the kernel which performs the
676  // triangular solve can assume a unit diagonal, take a short-cut
677  // and perform faster.
678 
679  L_->insertLocalValues (row_i, tmp_idx (), tmpv ());
680 #ifdef IFPACK2_WRITE_FACTORS
681  for (size_type ii = 0; ii < tmp_idx.size (); ++ii) {
682  ofsL << row_i << " " << tmp_idx[ii] << " " << tmpv[ii] << std::endl;
683  }
684 #endif
685 
686  tmp_idx.clear();
687  tmpv.clear();
688 
689  // Pick out the diagonal element, store its reciprocal.
690  if (cur_row[row_i] == zero) {
691  std::cerr << "Ifpack2::ILUT::Compute: zero pivot encountered! Replacing with rownorm and continuing...(You may need to set the parameter 'fact: absolute threshold'.)" << std::endl;
692  cur_row[row_i] = rownorm;
693  }
694  InvDiagU[row_i] = one / cur_row[row_i];
695 
696  // Non-inverted diagonal is stored for U:
697  tmp_idx.push_back(row_i);
698  tmpv.push_back(cur_row[row_i]);
699  unorm[row_i] = scalar_mag(cur_row[row_i]);
700  pattern[row_i] = UNUSED;
701 
702  // Now put indices and values for U into arrays and then into the U_ matrix.
703  // The first entry in U_cols is the diagonal, which we just handled, so we'll
704  // start our loop at j=1.
705 
706  size_type U_vals_heaplen = 0;
707  for(size_type j=1; j<U_cols.size(); ++j) {
708  local_ordinal_type col = U_cols[j];
709  if (pattern[col] != ORIG) {
710  if (U_vals_heaplen < fillU) {
711  add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
712  }
713  else if (U_vals_heaplen!=0 && scalar_mag(cur_row[col]) >
714  scalar_mag(cur_row[U_vals_heap.front()])) {
715  rm_heap_root(U_vals_heap, U_vals_heaplen, vals_comp);
716  add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
717  }
718  }
719  else {
720  tmp_idx.push_back(col);
721  tmpv.push_back(cur_row[col]);
722  unorm[row_i] += scalar_mag(cur_row[col]);
723  }
724  pattern[col] = UNUSED;
725  }
726 
727  for(size_type j=0; j<U_vals_heaplen; ++j) {
728  tmp_idx.push_back(U_vals_heap[j]);
729  tmpv.push_back(cur_row[U_vals_heap[j]]);
730  unorm[row_i] += scalar_mag(cur_row[U_vals_heap[j]]);
731  }
732 
733  unorm[row_i] /= (orig_U_len + U_vals_heaplen);
734 
735  U_->insertLocalValues(row_i, tmp_idx(), tmpv() );
736 #ifdef IFPACK2_WRITE_FACTORS
737  for(int ii=0; ii<tmp_idx.size(); ++ii) {
738  ofsU <<row_i<< " " <<tmp_idx[ii]<< " " <<tmpv[ii]<< std::endl;
739  }
740 #endif
741  tmp_idx.clear();
742  tmpv.clear();
743 
744  U_->getLocalRowView(row_i, Uindices[row_i], Ucoefs[row_i] );
745 
746  L_cols_heap.clear();
747  U_cols.clear();
748  L_vals_heap.clear();
749  U_vals_heap.clear();
750  } // end of for(row_i) loop
751 
752  // FIXME (mfh 03 Apr 2013) Do we need to supply a domain and range Map?
753  L_->fillComplete();
754  U_->fillComplete();
755  }
756  ComputeTime_ += timer.totalElapsedTime ();
757  IsComputed_ = true;
758  ++NumCompute_;
759 }
760 
761 
762 template <class MatrixType>
764 apply (const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
765  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
766  Teuchos::ETransp mode,
767  scalar_type alpha,
768  scalar_type beta) const
769 {
770  using Teuchos::RCP;
771  using Teuchos::rcp;
772  using Teuchos::rcpFromRef;
773  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
774 
775  Teuchos::Time timer ("ILUT::apply");
776  { // Timer scope for timing apply()
777  Teuchos::TimeMonitor timeMon (timer, true);
778 
779  TEUCHOS_TEST_FOR_EXCEPTION(
780  ! isComputed (), std::runtime_error,
781  "Ifpack2::ILUT::apply: You must call compute() to compute the incomplete "
782  "factorization, before calling apply().");
783 
784  TEUCHOS_TEST_FOR_EXCEPTION(
785  X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
786  "Ifpack2::ILUT::apply: X and Y must have the same number of columns. "
787  "X has " << X.getNumVectors () << " columns, but Y has "
788  << Y.getNumVectors () << " columns.");
789 
790  if (alpha == Teuchos::ScalarTraits<scalar_type>::zero ()) {
791  // alpha == 0, so we don't need to apply the operator.
792  //
793  // The special case for beta == 0 ensures that if Y contains Inf
794  // or NaN values, we replace them with 0 (following BLAS
795  // convention), rather than multiplying them by 0 to get NaN.
796  if (beta == Teuchos::ScalarTraits<scalar_type>::zero ()) {
797  Y.putScalar (beta);
798  } else {
799  Y.scale (beta);
800  }
801  return;
802  }
803 
804  // If beta != 0, create a temporary multivector Y_temp to hold the
805  // contents of alpha*M^{-1}*X. Otherwise, alias Y_temp to Y.
806  RCP<MV> Y_temp;
807  if (beta == Teuchos::ScalarTraits<scalar_type>::zero ()) {
808  Y_temp = rcpFromRef (Y);
809  } else {
810  Y_temp = rcp (new MV (Y.getMap (), Y.getNumVectors ()));
811  }
812 
813  // If X and Y are pointing to the same memory location, create an
814  // auxiliary vector, X_temp, so that we don't clobber the input
815  // when computing the output. Otherwise, alias X_temp to X.
816  RCP<const MV> X_temp;
817  {
818  auto X_lcl_host = X.template getLocalView<Kokkos::HostSpace> ();
819  auto Y_lcl_host = Y.template getLocalView<Kokkos::HostSpace> ();
820  if (X_lcl_host.ptr_on_device () == Y_lcl_host.ptr_on_device ()) {
821  X_temp = rcp (new MV (X, Teuchos::Copy));
822  } else {
823  X_temp = rcpFromRef (X);
824  }
825  }
826 
827  // Create a temporary multivector Y_mid to hold the intermediate
828  // between the L and U (or U and L, for the transpose or conjugate
829  // transpose case) solves.
830  RCP<MV> Y_mid = rcp (new MV (Y.getMap (), Y.getNumVectors ()));
831 
832  if (mode == Teuchos::NO_TRANS) { // Solve L U Y = X
833  L_->template localSolve<scalar_type, scalar_type> (*X_temp, *Y_mid, mode);
834  // FIXME (mfh 20 Aug 2013) Is it OK to use Y_temp for both the
835  // input and the output?
836  U_->template localSolve<scalar_type, scalar_type> (*Y_mid, *Y_temp, mode);
837  }
838  else { // Solve U^* L^* Y = X
839  U_->template localSolve<scalar_type, scalar_type> (*X_temp, *Y_mid, mode);
840  // FIXME (mfh 20 Aug 2013) Is it OK to use Y_temp for both the
841  // input and the output?
842  L_->template localSolve<scalar_type, scalar_type> (*Y_mid, *Y_temp, mode);
843  }
844 
845  if (beta == Teuchos::ScalarTraits<scalar_type>::zero ()) {
846  Y.scale (alpha);
847  } else { // beta != 0
848  Y.update (alpha, *Y_temp, beta);
849  }
850  }
851  ++NumApply_;
852  ApplyTime_ += timer.totalElapsedTime ();
853 }
854 
855 
856 template <class MatrixType>
857 std::string ILUT<MatrixType>::description () const
858 {
859  std::ostringstream os;
860 
861  // Output is a valid YAML dictionary in flow style. If you don't
862  // like everything on a single line, you should call describe()
863  // instead.
864  os << "\"Ifpack2::ILUT\": {";
865  os << "Initialized: " << (isInitialized () ? "true" : "false") << ", "
866  << "Computed: " << (isComputed () ? "true" : "false") << ", ";
867 
868  os << "Level-of-fill: " << getLevelOfFill() << ", "
869  << "absolute threshold: " << getAbsoluteThreshold() << ", "
870  << "relative threshold: " << getRelativeThreshold() << ", "
871  << "relaxation value: " << getRelaxValue() << ", ";
872 
873  if (A_.is_null ()) {
874  os << "Matrix: null";
875  }
876  else {
877  os << "Global matrix dimensions: ["
878  << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"
879  << ", Global nnz: " << A_->getGlobalNumEntries();
880  }
881 
882  os << "}";
883  return os.str ();
884 }
885 
886 
887 template <class MatrixType>
888 void
890 describe (Teuchos::FancyOStream& out,
891  const Teuchos::EVerbosityLevel verbLevel) const
892 {
893  using Teuchos::Comm;
894  using Teuchos::OSTab;
895  using Teuchos::RCP;
896  using Teuchos::TypeNameTraits;
897  using std::endl;
898  using Teuchos::VERB_DEFAULT;
899  using Teuchos::VERB_NONE;
900  using Teuchos::VERB_LOW;
901  using Teuchos::VERB_MEDIUM;
902  using Teuchos::VERB_HIGH;
903  using Teuchos::VERB_EXTREME;
904 
905  const Teuchos::EVerbosityLevel vl =
906  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
907  OSTab tab0 (out);
908 
909  if (vl > VERB_NONE) {
910  out << "\"Ifpack2::ILUT\":" << endl;
911  OSTab tab1 (out);
912  out << "MatrixType: " << TypeNameTraits<MatrixType>::name () << endl;
913  if (this->getObjectLabel () != "") {
914  out << "Label: \"" << this->getObjectLabel () << "\"" << endl;
915  }
916  out << "Initialized: " << (isInitialized () ? "true" : "false")
917  << endl
918  << "Computed: " << (isComputed () ? "true" : "false")
919  << endl
920  << "Level of fill: " << getLevelOfFill () << endl
921  << "Absolute threshold: " << getAbsoluteThreshold () << endl
922  << "Relative threshold: " << getRelativeThreshold () << endl
923  << "Relax value: " << getRelaxValue () << endl;
924 
925  if (isComputed () && vl >= VERB_HIGH) {
926  const double fillFraction =
927  (double) getGlobalNumEntries () / (double) A_->getGlobalNumEntries ();
928  const double nnzToRows =
929  (double) getGlobalNumEntries () / (double) U_->getGlobalNumRows ();
930 
931  out << "Dimensions of L: [" << L_->getGlobalNumRows () << ", "
932  << L_->getGlobalNumRows () << "]" << endl
933  << "Dimensions of U: [" << U_->getGlobalNumRows () << ", "
934  << U_->getGlobalNumRows () << "]" << endl
935  << "Number of nonzeros in factors: " << getGlobalNumEntries () << endl
936  << "Fill fraction of factors over A: " << fillFraction << endl
937  << "Ratio of nonzeros to rows: " << nnzToRows << endl;
938  }
939 
940  out << "Number of initialize calls: " << getNumInitialize () << endl
941  << "Number of compute calls: " << getNumCompute () << endl
942  << "Number of apply calls: " << getNumApply () << endl
943  << "Total time in seconds for initialize: " << getInitializeTime () << endl
944  << "Total time in seconds for compute: " << getComputeTime () << endl
945  << "Total time in seconds for apply: " << getApplyTime () << endl;
946 
947  out << "Local matrix:" << endl;
948  A_local_->describe (out, vl);
949  }
950 }
951 
952 template <class MatrixType>
953 Teuchos::RCP<const typename ILUT<MatrixType>::row_matrix_type>
954 ILUT<MatrixType>::makeLocalFilter (const Teuchos::RCP<const row_matrix_type>& A)
955 {
956  if (A->getComm ()->getSize () > 1) {
957  return Teuchos::rcp (new LocalFilter<row_matrix_type> (A));
958  } else {
959  return A;
960  }
961 }
962 
963 }//namespace Ifpack2
964 
965 
966 // FIXME (mfh 16 Sep 2014) We should really only use RowMatrix here!
967 // There's no need to instantiate for CrsMatrix too. All Ifpack2
968 // preconditioners can and should do dynamic casts if they need a type
969 // more specific than RowMatrix.
970 
971 #define IFPACK2_ILUT_INSTANT(S,LO,GO,N) \
972  template class Ifpack2::ILUT< Tpetra::RowMatrix<S, LO, GO, N> >;
973 
974 #endif /* IFPACK2_ILUT_DEF_HPP */
magnitude_type getDropTolerance() const
Gets the dropping tolerance.
Definition: Ifpack2_ILUT_decl.hpp:333
ILUT(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_ILUT_def.hpp:116
bool hasTransposeApply() const
Whether this object&#39;s apply() method can apply the transpose (or conjugate transpose, if applicable).
Definition: Ifpack2_ILUT_def.hpp:308
virtual ~ILUT()
Destructor.
Definition: Ifpack2_ILUT_def.hpp:134
global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_ILUT_def.hpp:350
void initialize()
Clear any previously computed factors.
Definition: Ifpack2_ILUT_def.hpp:389
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_ILUT_def.hpp:890
std::string description() const
Return a simple one-line description of this object.
Definition: Ifpack2_ILUT_def.hpp:857
ILUT (incomplete LU factorization with threshold) of a Tpetra sparse matrix.
Definition: Ifpack2_ILUT_decl.hpp:91
size_t getNodeNumEntries() const
Returns the number of nonzero entries in the local graph.
Definition: Ifpack2_ILUT_def.hpp:356
magnitude_type getRelaxValue() const
Get the relax value.
Definition: Ifpack2_ILUT_decl.hpp:328
void rm_heap_root(Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition: Ifpack2_Heap.hpp:92
Teuchos::RCP< const map_type > getRangeMap() const
Tpetra::Map representing the range of this operator.
Definition: Ifpack2_ILUT_def.hpp:297
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition: Ifpack2_ILUT_decl.hpp:323
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_ILUT_decl.hpp:106
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack2_ILUT_decl.hpp:318
void compute()
Compute factors L and U using the specified diagonal perturbation thresholds and relaxation parameter...
Definition: Ifpack2_ILUT_def.hpp:426
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_ILUT_def.hpp:362
int getNumInitialize() const
Returns the number of calls to Initialize().
Definition: Ifpack2_ILUT_def.hpp:314
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_ILUT_def.hpp:344
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the input matrix&#39;s communicator.
Definition: Ifpack2_ILUT_def.hpp:267
Definition: Ifpack2_Details_Amesos2Wrapper_decl.hpp:63
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Type of the Tpetra::CrsMatrix specialization that this class uses for the L and U factors...
Definition: Ifpack2_ILUT_decl.hpp:126
Teuchos::RCP< const map_type > getDomainMap() const
Tpetra::Map representing the domain of this operator.
Definition: Ifpack2_ILUT_def.hpp:285
int getLevelOfFill() const
The level of fill.
Definition: Ifpack2_ILUT_decl.hpp:313
bool isComputed() const
If compute() is completed, this query returns true, otherwise it returns false.
Definition: Ifpack2_ILUT_decl.hpp:215
void add_to_heap(const Ordinal &idx, Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition: Ifpack2_Heap.hpp:70
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:160
double getInitializeTime() const
Returns the time spent in Initialize().
Definition: Ifpack2_ILUT_def.hpp:332
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_ILUT_decl.hpp:200
int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_ILUT_def.hpp:326
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the ILUT preconditioner to X, resulting in Y.
Definition: Ifpack2_ILUT_def.hpp:764
double getComputeTime() const
Returns the time spent in Compute().
Definition: Ifpack2_ILUT_def.hpp:338
void setParameters(const Teuchos::ParameterList &params)
Set preconditioner parameters.
Definition: Ifpack2_ILUT_def.hpp:138
int getNumCompute() const
Returns the number of calls to Compute().
Definition: Ifpack2_ILUT_def.hpp:320
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_ILUT_decl.hpp:109
Teuchos::RCP< const row_matrix_type > getMatrix() const
Returns a reference to the matrix to be preconditioned.
Definition: Ifpack2_ILUT_def.hpp:278