Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_ILUT_def.hpp
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK2_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_LocalSparseTriangularSolver.hpp"
54 #include "Ifpack2_Parameters.hpp"
55 #include "Tpetra_CrsMatrix.hpp"
56 #include "Teuchos_Time.hpp"
58 #include <type_traits>
59 
60 namespace Ifpack2 {
61 
62  namespace {
63 
88  template<class ScalarType>
90  ilutDefaultDropTolerance () {
91  using Teuchos::as;
93  typedef typename STS::magnitudeType magnitude_type;
95 
96  // 1/2. Hopefully this can be represented in magnitude_type.
97  const magnitude_type oneHalf = STM::one() / (STM::one() + STM::one());
98 
99  // The min ensures that in case magnitude_type has very low
100  // precision, we'll at least get some value strictly less than
101  // one.
102  return std::min (as<magnitude_type> (1000) * STS::magnitude (STS::eps ()), oneHalf);
103  }
104 
105  // Full specialization for ScalarType = double.
106  // This specialization preserves ILUT's previous default behavior.
107  template<>
109  ilutDefaultDropTolerance<double> () {
110  return 1e-12;
111  }
112 
113  } // namespace (anonymous)
114 
115 
116 template <class MatrixType>
118  A_ (A),
119  Athresh_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
120  Rthresh_ (Teuchos::ScalarTraits<magnitude_type>::one ()),
121  RelaxValue_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
122  LevelOfFill_ (1.0),
123  DropTolerance_ (ilutDefaultDropTolerance<scalar_type> ()),
124  InitializeTime_ (0.0),
125  ComputeTime_ (0.0),
126  ApplyTime_ (0.0),
127  NumInitialize_ (0),
128  NumCompute_ (0),
129  NumApply_ (0),
130  IsInitialized_ (false),
131  IsComputed_ (false)
132 {
133  allocateSolvers();
134 }
135 
136 template <class MatrixType>
138 {}
139 
140 template<class MatrixType>
142 {
145 }
146 
147 namespace { // (anonymous)
148 
149  template<class MagnitudeType>
150  std::pair<MagnitudeType, bool>
151  getMagnitudeOrDoubleParameterAsMagnitude (const Teuchos::ParameterList& params,
152  const char parameterName[],
153  const MagnitudeType& currentValue)
154  {
155  const std::string pname (parameterName);
156 
157  if (params.isType<MagnitudeType> (pname)) {
158  const MagnitudeType value = params.get<MagnitudeType> (pname);
159  return {value, true};
160  }
161  else if (! std::is_same<MagnitudeType, double>::value &&
162  params.isType<double> (pname)) {
163  const MagnitudeType value = params.get<double> (pname);
164  return {value, true};
165  }
166  else {
167  return {currentValue, false};
168  }
169  }
170 
171 } // namespace (anonymous)
172 
173 
174 template <class MatrixType>
176 {
177  // Don't actually change the instance variables until we've checked
178  // all parameters. This ensures that setParameters satisfies the
179  // strong exception guarantee (i.e., is transactional).
180 
181  // Fill level in ILUT is a double, not a magnitude_type, because it
182  // depends on LO and GO, not on Scalar.
183  double fillLevel = LevelOfFill_;
184  bool gotFillLevel = false;
185  {
186  const std::string fillLevelParamName ("fact: ilut level-of-fill");
187  if (params.isType<double> (fillLevelParamName)) {
188  fillLevel = params.get<double> (fillLevelParamName);
189  gotFillLevel = true;
190  }
191  else if (! std::is_same<magnitude_type, double>::value &&
192  params.isType<double> (fillLevelParamName)) {
193  fillLevel = params.get<double> (fillLevelParamName);
194  gotFillLevel = true;
195  }
197  (gotFillLevel && fillLevel < 1.0, std::runtime_error,
198  "Ifpack2::ILUT: The \"fact: ilut level-of-fill\" parameter must be >= "
199  "1.0, but you set it to " << fillLevel << ". For ILUT, the fill level "
200  "means something different than it does for ILU(k). ILU(0) produces "
201  "factors with the same sparsity structure as the input matrix A. For "
202  "ILUT, level-of-fill = 1.0 will produce factors with nonzeros matching "
203  "the sparsity structure of A. level-of-fill > 1.0 allows for additional "
204  "fill-in.");
205  }
206 
207  const auto absThreshResult =
208  getMagnitudeOrDoubleParameterAsMagnitude (params, "fact: absolute threshold", Athresh_);
209  const auto relThreshResult =
210  getMagnitudeOrDoubleParameterAsMagnitude (params, "fact: relative threshold", Rthresh_);
211  const auto relaxResult =
212  getMagnitudeOrDoubleParameterAsMagnitude (params, "fact: relax value", RelaxValue_);
213  const auto dropResult =
214  getMagnitudeOrDoubleParameterAsMagnitude (params, "fact: drop tolerance", DropTolerance_);
215 
216  // Forward to trisolvers.
217  L_solver_->setParameters(params);
218  U_solver_->setParameters(params);
219 
220  if (gotFillLevel) {
221  LevelOfFill_ = fillLevel;
222  }
223  if (absThreshResult.second) {
224  Athresh_ = absThreshResult.first;
225  }
226  if (relThreshResult.second) {
227  Rthresh_ = relThreshResult.first;
228  }
229  if (relaxResult.second) {
230  RelaxValue_ = relaxResult.first;
231  }
232  if (dropResult.second) {
233  DropTolerance_ = dropResult.first;
234  }
235 }
236 
237 
238 template <class MatrixType>
242  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getComm: "
243  "The matrix is null. Please call setMatrix() with a nonnull input "
244  "before calling this method.");
245  return A_->getComm ();
246 }
247 
248 
249 template <class MatrixType>
252  return A_;
253 }
254 
255 
256 template <class MatrixType>
259 {
261  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getDomainMap: "
262  "The matrix is null. Please call setMatrix() with a nonnull input "
263  "before calling this method.");
264  return A_->getDomainMap ();
265 }
266 
267 
268 template <class MatrixType>
271 {
273  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getRangeMap: "
274  "The matrix is null. Please call setMatrix() with a nonnull input "
275  "before calling this method.");
276  return A_->getRangeMap ();
277 }
278 
279 
280 template <class MatrixType>
282  return true;
283 }
284 
285 
286 template <class MatrixType>
288  return NumInitialize_;
289 }
290 
291 
292 template <class MatrixType>
294  return NumCompute_;
295 }
296 
297 
298 template <class MatrixType>
300  return NumApply_;
301 }
302 
303 
304 template <class MatrixType>
306  return InitializeTime_;
307 }
308 
309 
310 template<class MatrixType>
312  return ComputeTime_;
313 }
314 
315 
316 template<class MatrixType>
318  return ApplyTime_;
319 }
320 
321 
322 template<class MatrixType>
325  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getNodeSmootherComplexity: "
326  "The input matrix A is null. Please call setMatrix() with a nonnull "
327  "input matrix, then call compute(), before calling this method.");
328  // ILUT methods cost roughly one apply + the nnz in the upper+lower triangles
329  return A_->getNodeNumEntries() + getNodeNumEntries();
330 }
331 
332 
333 template<class MatrixType>
335  return L_->getGlobalNumEntries () + U_->getGlobalNumEntries ();
336 }
337 
338 
339 template<class MatrixType>
341  return L_->getNodeNumEntries () + U_->getNodeNumEntries ();
342 }
343 
344 
345 template<class MatrixType>
347 {
348  if (A.getRawPtr () != A_.getRawPtr ()) {
349  // Check in serial or one-process mode if the matrix is square.
351  ! A.is_null () && A->getComm ()->getSize () == 1 &&
352  A->getNodeNumRows () != A->getNodeNumCols (),
353  std::runtime_error, "Ifpack2::ILUT::setMatrix: If A's communicator only "
354  "contains one process, then A must be square. Instead, you provided a "
355  "matrix A with " << A->getNodeNumRows () << " rows and "
356  << A->getNodeNumCols () << " columns.");
357 
358  // It's legal for A to be null; in that case, you may not call
359  // initialize() until calling setMatrix() with a nonnull input.
360  // Regardless, setting the matrix invalidates any previous
361  // factorization.
362  IsInitialized_ = false;
363  IsComputed_ = false;
364  A_local_ = Teuchos::null;
365 
366  // The sparse triangular solvers get a triangular factor as their
367  // input matrix. The triangular factors L_ and U_ are getting
368  // reset, so we reset the solvers' matrices to null. Do that
369  // before setting L_ and U_ to null, so that latter step actually
370  // frees the factors.
371  if (! L_solver_.is_null ()) {
372  L_solver_->setMatrix (Teuchos::null);
373  }
374  if (! U_solver_.is_null ()) {
375  U_solver_->setMatrix (Teuchos::null);
376  }
377 
378  L_ = Teuchos::null;
379  U_ = Teuchos::null;
380  A_ = A;
381  }
382 }
383 
384 
385 template<class MatrixType>
387 {
388  Teuchos::Time timer ("ILUT::initialize");
389  {
390  Teuchos::TimeMonitor timeMon (timer);
391 
392  // Check that the matrix is nonnull.
394  A_.is_null (), std::runtime_error, "Ifpack2::ILUT::initialize: "
395  "The matrix to precondition is null. Please call setMatrix() with a "
396  "nonnull input before calling this method.");
397 
398  // Clear any previous computations.
399  IsInitialized_ = false;
400  IsComputed_ = false;
401  A_local_ = Teuchos::null;
402  L_ = Teuchos::null;
403  U_ = Teuchos::null;
404 
405  A_local_ = makeLocalFilter (A_); // Compute the local filter.
406 
407  IsInitialized_ = true;
408  ++NumInitialize_;
409  }
410  InitializeTime_ += timer.totalElapsedTime ();
411 }
412 
413 
414 template<typename ScalarType>
416 scalar_mag (const ScalarType& s)
417 {
419 }
420 
421 
422 template<class MatrixType>
424 {
425  using Teuchos::Array;
426  using Teuchos::ArrayRCP;
427  using Teuchos::ArrayView;
428  using Teuchos::as;
429  using Teuchos::rcp;
430  using Teuchos::reduceAll;
431 
432  //--------------------------------------------------------------------------
433  // Ifpack2::ILUT is a translation of the Aztec ILUT implementation. The Aztec
434  // ILUT implementation was written by Ray Tuminaro.
435  //
436  // This isn't an exact translation of the Aztec ILUT algorithm, for the
437  // following reasons:
438  // 1. Minor differences result from the fact that Aztec factors a MSR format
439  // matrix in place, while the code below factors an input CrsMatrix which
440  // remains untouched and stores the resulting factors in separate L and U
441  // CrsMatrix objects.
442  // Also, the Aztec code begins by shifting the matrix pointers back
443  // by one, and the pointer contents back by one, and then using 1-based
444  // Fortran-style indexing in the algorithm. This Ifpack2 code uses C-style
445  // 0-based indexing throughout.
446  // 2. Aztec stores the inverse of the diagonal of U. This Ifpack2 code
447  // stores the non-inverted diagonal in U.
448  // The triangular solves (in Ifpack2::ILUT::apply()) are performed by
449  // calling the Tpetra::CrsMatrix::solve method on the L and U objects, and
450  // this requires U to contain the non-inverted diagonal.
451  //
452  // ABW.
453  //--------------------------------------------------------------------------
454 
455  // Don't count initialization in the compute() time.
456  if (! isInitialized ()) {
457  initialize ();
458  }
459 
460  Teuchos::Time timer ("ILUT::compute");
461  { // Timer scope for timing compute()
462  Teuchos::TimeMonitor timeMon (timer, true);
463  const scalar_type zero = STS::zero ();
464  const scalar_type one = STS::one ();
465 
466  const local_ordinal_type myNumRows = A_local_->getNodeNumRows ();
467  L_ = rcp (new crs_matrix_type (A_local_->getRowMap (), A_local_->getColMap (), 0));
468  U_ = rcp (new crs_matrix_type (A_local_->getRowMap (), A_local_->getColMap (), 0));
469 
470  // CGB: note, this caching approach may not be necessary anymore
471  // We will store ArrayView objects that are views of the rows of U, so that
472  // we don't have to repeatedly retrieve the view for each row. These will
473  // be populated row by row as the factorization proceeds.
474  Array<ArrayView<const local_ordinal_type> > Uindices (myNumRows);
475  Array<ArrayView<const scalar_type> > Ucoefs (myNumRows);
476 
477  // If this macro is defined, files containing the L and U factors
478  // will be written. DON'T CHECK IN THE CODE WITH THIS MACRO ENABLED!!!
479  // #define IFPACK2_WRITE_FACTORS
480 #ifdef IFPACK2_WRITE_FACTORS
481  std::ofstream ofsL("L.tif.mtx", std::ios::out);
482  std::ofstream ofsU("U.tif.mtx", std::ios::out);
483 #endif
484 
485  // Calculate how much fill will be allowed in addition to the
486  // space that corresponds to the input matrix entries.
487  double local_nnz = static_cast<double> (A_local_->getNodeNumEntries ());
488  double fill = ((getLevelOfFill () - 1.0) * local_nnz) / (2 * myNumRows);
489 
490  // std::ceil gives the smallest integer larger than the argument.
491  // this may give a slightly different result than Aztec's fill value in
492  // some cases.
493  double fill_ceil=std::ceil(fill);
494 
495  // Similarly to Aztec, we will allow the same amount of fill for each
496  // row, half in L and half in U.
497  size_type fillL = static_cast<size_type>(fill_ceil);
498  size_type fillU = static_cast<size_type>(fill_ceil);
499 
500  Array<scalar_type> InvDiagU (myNumRows, zero);
501 
502  Array<local_ordinal_type> tmp_idx;
503  Array<scalar_type> tmpv;
504 
505  enum { UNUSED, ORIG, FILL };
506  local_ordinal_type max_col = myNumRows;
507 
508  Array<int> pattern(max_col, UNUSED);
509  Array<scalar_type> cur_row(max_col, zero);
510  Array<magnitude_type> unorm(max_col);
511  magnitude_type rownorm;
512  Array<local_ordinal_type> L_cols_heap;
513  Array<local_ordinal_type> U_cols;
514  Array<local_ordinal_type> L_vals_heap;
515  Array<local_ordinal_type> U_vals_heap;
516 
517  // A comparison object which will be used to create 'heaps' of indices
518  // that are ordered according to the corresponding values in the
519  // 'cur_row' array.
520  greater_indirect<scalar_type,local_ordinal_type> vals_comp(cur_row);
521 
522  // =================== //
523  // start factorization //
524  // =================== //
525 
526  ArrayRCP<local_ordinal_type> ColIndicesARCP;
527  ArrayRCP<scalar_type> ColValuesARCP;
528  if (! A_local_->supportsRowViews ()) {
529  const size_t maxnz = A_local_->getNodeMaxNumRowEntries ();
530  ColIndicesARCP.resize (maxnz);
531  ColValuesARCP.resize (maxnz);
532  }
533 
534  for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
535  ArrayView<const local_ordinal_type> ColIndicesA;
536  ArrayView<const scalar_type> ColValuesA;
537  size_t RowNnz;
538 
539  if (A_local_->supportsRowViews ()) {
540  A_local_->getLocalRowView (row_i, ColIndicesA, ColValuesA);
541  RowNnz = ColIndicesA.size ();
542  }
543  else {
544  A_local_->getLocalRowCopy (row_i, ColIndicesARCP (), ColValuesARCP (), RowNnz);
545  ColIndicesA = ColIndicesARCP (0, RowNnz);
546  ColValuesA = ColValuesARCP (0, RowNnz);
547  }
548 
549  // Always include the diagonal in the U factor. The value should get
550  // set in the next loop below.
551  U_cols.push_back(row_i);
552  cur_row[row_i] = zero;
553  pattern[row_i] = ORIG;
554 
555  size_type L_cols_heaplen = 0;
556  rownorm = STM::zero ();
557  for (size_t i = 0; i < RowNnz; ++i) {
558  if (ColIndicesA[i] < myNumRows) {
559  if (ColIndicesA[i] < row_i) {
560  add_to_heap(ColIndicesA[i], L_cols_heap, L_cols_heaplen);
561  }
562  else if (ColIndicesA[i] > row_i) {
563  U_cols.push_back(ColIndicesA[i]);
564  }
565 
566  cur_row[ColIndicesA[i]] = ColValuesA[i];
567  pattern[ColIndicesA[i]] = ORIG;
568  rownorm += scalar_mag(ColValuesA[i]);
569  }
570  }
571 
572  // Alter the diagonal according to the absolute-threshold and
573  // relative-threshold values. If not set, those values default
574  // to zero and one respectively.
575  const magnitude_type rthresh = getRelativeThreshold();
576  const scalar_type& v = cur_row[row_i];
577  cur_row[row_i] = as<scalar_type> (getAbsoluteThreshold() * IFPACK2_SGN(v)) + rthresh*v;
578 
579  size_type orig_U_len = U_cols.size();
580  RowNnz = L_cols_heap.size() + orig_U_len;
581  rownorm = getDropTolerance() * rownorm/RowNnz;
582 
583  // The following while loop corresponds to the 'L30' goto's in Aztec.
584  size_type L_vals_heaplen = 0;
585  while (L_cols_heaplen > 0) {
586  local_ordinal_type row_k = L_cols_heap.front();
587 
588  scalar_type multiplier = cur_row[row_k] * InvDiagU[row_k];
589  cur_row[row_k] = multiplier;
590  magnitude_type mag_mult = scalar_mag(multiplier);
591  if (mag_mult*unorm[row_k] < rownorm) {
592  pattern[row_k] = UNUSED;
593  rm_heap_root(L_cols_heap, L_cols_heaplen);
594  continue;
595  }
596  if (pattern[row_k] != ORIG) {
597  if (L_vals_heaplen < fillL) {
598  add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
599  }
600  else if (L_vals_heaplen==0 ||
601  mag_mult < scalar_mag(cur_row[L_vals_heap.front()])) {
602  pattern[row_k] = UNUSED;
603  rm_heap_root(L_cols_heap, L_cols_heaplen);
604  continue;
605  }
606  else {
607  pattern[L_vals_heap.front()] = UNUSED;
608  rm_heap_root(L_vals_heap, L_vals_heaplen, vals_comp);
609  add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
610  }
611  }
612 
613  /* Reduce current row */
614 
615  ArrayView<const local_ordinal_type>& ColIndicesU = Uindices[row_k];
616  ArrayView<const scalar_type>& ColValuesU = Ucoefs[row_k];
617  size_type ColNnzU = ColIndicesU.size();
618 
619  for(size_type j=0; j<ColNnzU; ++j) {
620  if (ColIndicesU[j] > row_k) {
621  scalar_type tmp = multiplier * ColValuesU[j];
622  local_ordinal_type col_j = ColIndicesU[j];
623  if (pattern[col_j] != UNUSED) {
624  cur_row[col_j] -= tmp;
625  }
626  else if (scalar_mag(tmp) > rownorm) {
627  cur_row[col_j] = -tmp;
628  pattern[col_j] = FILL;
629  if (col_j > row_i) {
630  U_cols.push_back(col_j);
631  }
632  else {
633  add_to_heap(col_j, L_cols_heap, L_cols_heaplen);
634  }
635  }
636  }
637  }
638 
639  rm_heap_root(L_cols_heap, L_cols_heaplen);
640  }//end of while(L_cols_heaplen) loop
641 
642 
643  // Put indices and values for L into arrays and then into the L_ matrix.
644 
645  // first, the original entries from the L section of A:
646  for (size_type i = 0; i < ColIndicesA.size (); ++i) {
647  if (ColIndicesA[i] < row_i) {
648  tmp_idx.push_back(ColIndicesA[i]);
649  tmpv.push_back(cur_row[ColIndicesA[i]]);
650  pattern[ColIndicesA[i]] = UNUSED;
651  }
652  }
653 
654  // next, the L entries resulting from fill:
655  for (size_type j = 0; j < L_vals_heaplen; ++j) {
656  tmp_idx.push_back(L_vals_heap[j]);
657  tmpv.push_back(cur_row[L_vals_heap[j]]);
658  pattern[L_vals_heap[j]] = UNUSED;
659  }
660 
661  // L has a one on the diagonal, but we don't explicitly store
662  // it. If we don't store it, then the kernel which performs the
663  // triangular solve can assume a unit diagonal, take a short-cut
664  // and perform faster.
665 
666  L_->insertLocalValues (row_i, tmp_idx (), tmpv ());
667 #ifdef IFPACK2_WRITE_FACTORS
668  for (size_type ii = 0; ii < tmp_idx.size (); ++ii) {
669  ofsL << row_i << " " << tmp_idx[ii] << " " << tmpv[ii] << std::endl;
670  }
671 #endif
672 
673  tmp_idx.clear();
674  tmpv.clear();
675 
676  // Pick out the diagonal element, store its reciprocal.
677  if (cur_row[row_i] == zero) {
678  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;
679  cur_row[row_i] = rownorm;
680  }
681  InvDiagU[row_i] = one / cur_row[row_i];
682 
683  // Non-inverted diagonal is stored for U:
684  tmp_idx.push_back(row_i);
685  tmpv.push_back(cur_row[row_i]);
686  unorm[row_i] = scalar_mag(cur_row[row_i]);
687  pattern[row_i] = UNUSED;
688 
689  // Now put indices and values for U into arrays and then into the U_ matrix.
690  // The first entry in U_cols is the diagonal, which we just handled, so we'll
691  // start our loop at j=1.
692 
693  size_type U_vals_heaplen = 0;
694  for(size_type j=1; j<U_cols.size(); ++j) {
695  local_ordinal_type col = U_cols[j];
696  if (pattern[col] != ORIG) {
697  if (U_vals_heaplen < fillU) {
698  add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
699  }
700  else if (U_vals_heaplen!=0 && scalar_mag(cur_row[col]) >
701  scalar_mag(cur_row[U_vals_heap.front()])) {
702  rm_heap_root(U_vals_heap, U_vals_heaplen, vals_comp);
703  add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
704  }
705  }
706  else {
707  tmp_idx.push_back(col);
708  tmpv.push_back(cur_row[col]);
709  unorm[row_i] += scalar_mag(cur_row[col]);
710  }
711  pattern[col] = UNUSED;
712  }
713 
714  for(size_type j=0; j<U_vals_heaplen; ++j) {
715  tmp_idx.push_back(U_vals_heap[j]);
716  tmpv.push_back(cur_row[U_vals_heap[j]]);
717  unorm[row_i] += scalar_mag(cur_row[U_vals_heap[j]]);
718  }
719 
720  unorm[row_i] /= (orig_U_len + U_vals_heaplen);
721 
722  U_->insertLocalValues(row_i, tmp_idx(), tmpv() );
723 #ifdef IFPACK2_WRITE_FACTORS
724  for(int ii=0; ii<tmp_idx.size(); ++ii) {
725  ofsU <<row_i<< " " <<tmp_idx[ii]<< " " <<tmpv[ii]<< std::endl;
726  }
727 #endif
728  tmp_idx.clear();
729  tmpv.clear();
730 
731  U_->getLocalRowView(row_i, Uindices[row_i], Ucoefs[row_i] );
732 
733  L_cols_heap.clear();
734  U_cols.clear();
735  L_vals_heap.clear();
736  U_vals_heap.clear();
737  } // end of for(row_i) loop
738 
739  // FIXME (mfh 03 Apr 2013) Do we need to supply a domain and range Map?
740  L_->fillComplete();
741  U_->fillComplete();
742 
743  L_solver_->setMatrix(L_);
744  L_solver_->initialize ();
745  L_solver_->compute ();
746 
747  U_solver_->setMatrix(U_);
748  U_solver_->initialize ();
749  U_solver_->compute ();
750  }
751  ComputeTime_ += timer.totalElapsedTime ();
752  IsComputed_ = true;
753  ++NumCompute_;
754 }
755 
756 
757 template <class MatrixType>
759 apply (const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
760  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
761  Teuchos::ETransp mode,
762  scalar_type alpha,
763  scalar_type beta) const
764 {
765  using Teuchos::RCP;
766  using Teuchos::rcp;
767  using Teuchos::rcpFromRef;
769 
771  ! isComputed (), std::runtime_error,
772  "Ifpack2::ILUT::apply: You must call compute() to compute the incomplete "
773  "factorization, before calling apply().");
774 
776  X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
777  "Ifpack2::ILUT::apply: X and Y must have the same number of columns. "
778  "X has " << X.getNumVectors () << " columns, but Y has "
779  << Y.getNumVectors () << " columns.");
780 
781  const scalar_type one = STS::one ();
782  const scalar_type zero = STS::zero ();
783 
784  Teuchos::Time timer ("ILUT::apply");
785  { // Start timing
786  Teuchos::TimeMonitor timeMon (timer, true);
787 
788  if (alpha == one && beta == zero) {
789  if (mode == Teuchos::NO_TRANS) { // Solve L (U Y) = X for Y.
790  // Start by solving L Y = X for Y.
791  L_solver_->apply (X, Y, mode);
792 
793  // Solve U Y = Y.
794  U_solver_->apply (Y, Y, mode);
795  }
796  else { // Solve U^P (L^P Y)) = X for Y (where P is * or T).
797 
798  // Start by solving U^P Y = X for Y.
799  U_solver_->apply (X, Y, mode);
800 
801  // Solve L^P Y = Y.
802  L_solver_->apply (Y, Y, mode);
803  }
804  }
805  else { // alpha != 1 or beta != 0
806  if (alpha == zero) {
807  // The special case for beta == 0 ensures that if Y contains Inf
808  // or NaN values, we replace them with 0 (following BLAS
809  // convention), rather than multiplying them by 0 to get NaN.
810  if (beta == zero) {
811  Y.putScalar (zero);
812  } else {
813  Y.scale (beta);
814  }
815  } else { // alpha != zero
816  MV Y_tmp (Y.getMap (), Y.getNumVectors ());
817  apply (X, Y_tmp, mode);
818  Y.update (alpha, Y_tmp, beta);
819  }
820  }
821  }//end timing
822 
823  ++NumApply_;
824  ApplyTime_ += timer.totalElapsedTime ();
825 }
826 
827 
828 template <class MatrixType>
829 std::string ILUT<MatrixType>::description () const
830 {
831  std::ostringstream os;
832 
833  // Output is a valid YAML dictionary in flow style. If you don't
834  // like everything on a single line, you should call describe()
835  // instead.
836  os << "\"Ifpack2::ILUT\": {";
837  os << "Initialized: " << (isInitialized () ? "true" : "false") << ", "
838  << "Computed: " << (isComputed () ? "true" : "false") << ", ";
839 
840  os << "Level-of-fill: " << getLevelOfFill() << ", "
841  << "absolute threshold: " << getAbsoluteThreshold() << ", "
842  << "relative threshold: " << getRelativeThreshold() << ", "
843  << "relaxation value: " << getRelaxValue() << ", ";
844 
845  if (A_.is_null ()) {
846  os << "Matrix: null";
847  }
848  else {
849  os << "Global matrix dimensions: ["
850  << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"
851  << ", Global nnz: " << A_->getGlobalNumEntries();
852  }
853 
854  os << "}";
855  return os.str ();
856 }
857 
858 
859 template <class MatrixType>
860 void
863  const Teuchos::EVerbosityLevel verbLevel) const
864 {
865  using Teuchos::Comm;
866  using Teuchos::OSTab;
867  using Teuchos::RCP;
869  using std::endl;
870  using Teuchos::VERB_DEFAULT;
871  using Teuchos::VERB_NONE;
872  using Teuchos::VERB_LOW;
873  using Teuchos::VERB_MEDIUM;
874  using Teuchos::VERB_HIGH;
875  using Teuchos::VERB_EXTREME;
876 
877  const Teuchos::EVerbosityLevel vl =
878  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
879  OSTab tab0 (out);
880 
881  if (vl > VERB_NONE) {
882  out << "\"Ifpack2::ILUT\":" << endl;
883  OSTab tab1 (out);
884  out << "MatrixType: " << TypeNameTraits<MatrixType>::name () << endl;
885  if (this->getObjectLabel () != "") {
886  out << "Label: \"" << this->getObjectLabel () << "\"" << endl;
887  }
888  out << "Initialized: " << (isInitialized () ? "true" : "false")
889  << endl
890  << "Computed: " << (isComputed () ? "true" : "false")
891  << endl
892  << "Level of fill: " << getLevelOfFill () << endl
893  << "Absolute threshold: " << getAbsoluteThreshold () << endl
894  << "Relative threshold: " << getRelativeThreshold () << endl
895  << "Relax value: " << getRelaxValue () << endl;
896 
897  if (isComputed () && vl >= VERB_HIGH) {
898  const double fillFraction =
899  (double) getGlobalNumEntries () / (double) A_->getGlobalNumEntries ();
900  const double nnzToRows =
901  (double) getGlobalNumEntries () / (double) U_->getGlobalNumRows ();
902 
903  out << "Dimensions of L: [" << L_->getGlobalNumRows () << ", "
904  << L_->getGlobalNumRows () << "]" << endl
905  << "Dimensions of U: [" << U_->getGlobalNumRows () << ", "
906  << U_->getGlobalNumRows () << "]" << endl
907  << "Number of nonzeros in factors: " << getGlobalNumEntries () << endl
908  << "Fill fraction of factors over A: " << fillFraction << endl
909  << "Ratio of nonzeros to rows: " << nnzToRows << endl;
910  }
911 
912  out << "Number of initialize calls: " << getNumInitialize () << endl
913  << "Number of compute calls: " << getNumCompute () << endl
914  << "Number of apply calls: " << getNumApply () << endl
915  << "Total time in seconds for initialize: " << getInitializeTime () << endl
916  << "Total time in seconds for compute: " << getComputeTime () << endl
917  << "Total time in seconds for apply: " << getApplyTime () << endl;
918 
919  out << "Local matrix:" << endl;
920  A_local_->describe (out, vl);
921  }
922 }
923 
924 template <class MatrixType>
927 {
928  if (A->getComm ()->getSize () > 1) {
930  } else {
931  return A;
932  }
933 }
934 
935 }//namespace Ifpack2
936 
937 
938 // FIXME (mfh 16 Sep 2014) We should really only use RowMatrix here!
939 // There's no need to instantiate for CrsMatrix too. All Ifpack2
940 // preconditioners can and should do dynamic casts if they need a type
941 // more specific than RowMatrix.
942 
943 #define IFPACK2_ILUT_INSTANT(S,LO,GO,N) \
944  template class Ifpack2::ILUT< Tpetra::RowMatrix<S, LO, GO, N> >;
945 
946 #endif /* IFPACK2_ILUT_DEF_HPP */
ILUT(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_ILUT_def.hpp:117
bool hasTransposeApply() const
Whether this object&#39;s apply() method can apply the transpose (or conjugate transpose, if applicable).
Definition: Ifpack2_ILUT_def.hpp:281
virtual ~ILUT()
Destructor.
Definition: Ifpack2_ILUT_def.hpp:137
basic_OSTab< char > OSTab
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_ILUT_def.hpp:334
void initialize()
Clear any previously computed factors.
Definition: Ifpack2_ILUT_def.hpp:386
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:862
std::string description() const
Return a simple one-line description of this object.
Definition: Ifpack2_ILUT_def.hpp:829
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:340
void rm_heap_root(Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition: Ifpack2_Heap.hpp:92
T * getRawPtr() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const map_type > getRangeMap() const
Tpetra::Map representing the range of this operator.
Definition: Ifpack2_ILUT_def.hpp:270
&quot;Preconditioner&quot; that solves local sparse triangular systems.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:77
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_ILUT_decl.hpp:106
TEUCHOS_DEPRECATED void reduceAll(const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, Packet *globalReduct)
void compute()
Compute factors L and U using the specified diagonal perturbation thresholds and relaxation parameter...
Definition: Ifpack2_ILUT_def.hpp:423
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_ILUT_def.hpp:346
static magnitudeType magnitude(T a)
int getNumInitialize() const
Returns the number of calls to Initialize().
Definition: Ifpack2_ILUT_def.hpp:287
double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_ILUT_def.hpp:317
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the input matrix&#39;s communicator.
Definition: Ifpack2_ILUT_def.hpp:240
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
TypeTo as(const TypeFrom &t)
Teuchos::RCP< const map_type > getDomainMap() const
Tpetra::Map representing the domain of this operator.
Definition: Ifpack2_ILUT_def.hpp:258
bool isType(const std::string &name) const
double totalElapsedTime(bool readCurrentTime=false) const
void add_to_heap(const Ordinal &idx, Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition: Ifpack2_Heap.hpp:70
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_ILUT_def.hpp:323
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:305
int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_ILUT_def.hpp:299
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:759
double getComputeTime() const
Returns the time spent in Compute().
Definition: Ifpack2_ILUT_def.hpp:311
void setParameters(const Teuchos::ParameterList &params)
Set preconditioner parameters.
Definition: Ifpack2_ILUT_def.hpp:175
int getNumCompute() const
Returns the number of calls to Compute().
Definition: Ifpack2_ILUT_def.hpp:293
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_ILUT_decl.hpp:109
bool is_null() const
Teuchos::RCP< const row_matrix_type > getMatrix() const
Returns a reference to the matrix to be preconditioned.
Definition: Ifpack2_ILUT_def.hpp:251