Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_RILUK_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_CRSRILUK_DEF_HPP
44 #define IFPACK2_CRSRILUK_DEF_HPP
45 
46 #include "Ifpack2_LocalFilter.hpp"
47 #include "Tpetra_CrsMatrix.hpp"
48 #include "Ifpack2_LocalSparseTriangularSolver.hpp"
49 
50 namespace Ifpack2 {
51 
52 template<class MatrixType>
54  : A_ (Matrix_in),
55  LevelOfFill_ (0),
56  isAllocated_ (false),
57  isInitialized_ (false),
58  isComputed_ (false),
59  numInitialize_ (0),
60  numCompute_ (0),
61  numApply_ (0),
62  initializeTime_ (0.0),
63  computeTime_ (0.0),
64  applyTime_ (0.0),
65  RelaxValue_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
66  Athresh_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
67  Rthresh_ (Teuchos::ScalarTraits<magnitude_type>::one ())
68 {
69  allocateSolvers();
70 }
71 
72 
73 template<class MatrixType>
75  : A_ (Matrix_in),
76  LevelOfFill_ (0),
77  isAllocated_ (false),
78  isInitialized_ (false),
79  isComputed_ (false),
80  numInitialize_ (0),
81  numCompute_ (0),
82  numApply_ (0),
83  initializeTime_ (0.0),
84  computeTime_ (0.0),
85  applyTime_ (0.0),
86  RelaxValue_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
87  Athresh_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
88  Rthresh_ (Teuchos::ScalarTraits<magnitude_type>::one ())
89 {
90  allocateSolvers();
91 }
92 
93 
94 template<class MatrixType>
96 
97 template<class MatrixType>
99 {
102 }
103 
104 template<class MatrixType>
105 void
107 {
108  // It's legal for A to be null; in that case, you may not call
109  // initialize() until calling setMatrix() with a nonnull input.
110  // Regardless, setting the matrix invalidates any previous
111  // factorization.
112  if (A.getRawPtr () != A_.getRawPtr ()) {
113  isAllocated_ = false;
114  isInitialized_ = false;
115  isComputed_ = false;
116  A_local_ = Teuchos::null;
117  Graph_ = Teuchos::null;
118 
119  // The sparse triangular solvers get a triangular factor as their
120  // input matrix. The triangular factors L_ and U_ are getting
121  // reset, so we reset the solvers' matrices to null. Do that
122  // before setting L_ and U_ to null, so that latter step actually
123  // frees the factors.
124  if (! L_solver_.is_null ()) {
125  L_solver_->setMatrix (Teuchos::null);
126  }
127  if (! U_solver_.is_null ()) {
128  U_solver_->setMatrix (Teuchos::null);
129  }
130 
131  L_ = Teuchos::null;
132  U_ = Teuchos::null;
133  D_ = Teuchos::null;
134  A_ = A;
135  }
136 }
137 
138 
139 
140 template<class MatrixType>
143 {
145  L_.is_null (), std::runtime_error, "Ifpack2::RILUK::getL: The L factor "
146  "is null. Please call initialize() (and preferably also compute()) "
147  "before calling this method. If the input matrix has not yet been set, "
148  "you must first call setMatrix() with a nonnull input matrix before you "
149  "may call initialize() or compute().");
150  return *L_;
151 }
152 
153 
154 template<class MatrixType>
155 const Tpetra::Vector<typename RILUK<MatrixType>::scalar_type,
160 {
162  D_.is_null (), std::runtime_error, "Ifpack2::RILUK::getD: The D factor "
163  "(of diagonal entries) is null. Please call initialize() (and "
164  "preferably also compute()) before calling this method. If the input "
165  "matrix has not yet been set, you must first call setMatrix() with a "
166  "nonnull input matrix before you may call initialize() or compute().");
167  return *D_;
168 }
169 
170 
171 template<class MatrixType>
174 {
176  U_.is_null (), std::runtime_error, "Ifpack2::RILUK::getU: The U factor "
177  "is null. Please call initialize() (and preferably also compute()) "
178  "before calling this method. If the input matrix has not yet been set, "
179  "you must first call setMatrix() with a nonnull input matrix before you "
180  "may call initialize() or compute().");
181  return *U_;
182 }
183 
184 
185 template<class MatrixType>
188  A_.is_null (), std::runtime_error, "Ifpack2::RILUK::getNodeSmootherComplexity: "
189  "The input matrix A is null. Please call setMatrix() with a nonnull "
190  "input matrix, then call compute(), before calling this method.");
191  // RILUK methods cost roughly one apply + the nnz in the upper+lower triangles
192  if(!L_.is_null() && !U_.is_null())
193  return A_->getNodeNumEntries() + L_->getNodeNumEntries() + U_->getNodeNumEntries();
194  else
195  return 0;
196 }
197 
198 
199 
200 template<class MatrixType>
203  typename RILUK<MatrixType>::node_type> >
206  A_.is_null (), std::runtime_error, "Ifpack2::RILUK::getDomainMap: "
207  "The matrix is null. Please call setMatrix() with a nonnull input "
208  "before calling this method.");
209 
210  // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
212  Graph_.is_null (), std::runtime_error, "Ifpack2::RILUK::getDomainMap: "
213  "The computed graph is null. Please call initialize() before calling "
214  "this method.");
215  return Graph_->getL_Graph ()->getDomainMap ();
216 }
217 
218 
219 template<class MatrixType>
222  typename RILUK<MatrixType>::node_type> >
225  A_.is_null (), std::runtime_error, "Ifpack2::RILUK::getRangeMap: "
226  "The matrix is null. Please call setMatrix() with a nonnull input "
227  "before calling this method.");
228 
229  // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
231  Graph_.is_null (), std::runtime_error, "Ifpack2::RILUK::getRangeMap: "
232  "The computed graph is null. Please call initialize() before calling "
233  "this method.");
234  return Graph_->getL_Graph ()->getRangeMap ();
235 }
236 
237 
238 template<class MatrixType>
240 {
241  using Teuchos::null;
242  using Teuchos::rcp;
243 
244  if (! isAllocated_) {
245  // Deallocate any existing storage. This avoids storing 2x
246  // memory, since RCP op= does not deallocate until after the
247  // assignment.
248  L_ = null;
249  U_ = null;
250  D_ = null;
251 
252  // Allocate Matrix using ILUK graphs
253  L_ = rcp (new crs_matrix_type (Graph_->getL_Graph ()));
254  U_ = rcp (new crs_matrix_type (Graph_->getU_Graph ()));
255  L_->setAllToScalar (STS::zero ()); // Zero out L and U matrices
256  U_->setAllToScalar (STS::zero ());
257 
258  // FIXME (mfh 24 Jan 2014) This assumes domain == range Map for L and U.
259  L_->fillComplete ();
260  U_->fillComplete ();
261  D_ = rcp (new vec_type (Graph_->getL_Graph ()->getRowMap ()));
262  }
263  isAllocated_ = true;
264 }
265 
266 
267 template<class MatrixType>
268 void
271 {
272  using Teuchos::as;
275 
276  // Default values of the various parameters.
277  int fillLevel = 0;
278  magnitude_type absThresh = STM::zero ();
279  magnitude_type relThresh = STM::one ();
280  magnitude_type relaxValue = STM::zero ();
281 
282  //
283  // "fact: iluk level-of-fill" parsing is more complicated, because
284  // we want to allow as many types as make sense. int is the native
285  // type, but we also want to accept magnitude_type (for
286  // compatibility with ILUT) and double (for backwards compatibilty
287  // with ILUT).
288  //
289 
290  bool gotFillLevel = false;
291  try {
292  fillLevel = params.get<int> ("fact: iluk level-of-fill");
293  gotFillLevel = true;
294  }
295  catch (InvalidParameterType&) {
296  // Throwing again in the catch block would just unwind the stack.
297  // Instead, we do nothing here, and check the Boolean outside to
298  // see if we got the value.
299  }
300  catch (InvalidParameterName&) {
301  gotFillLevel = true; // Accept the default value.
302  }
303 
304  if (! gotFillLevel) {
305  try {
306  // Try global_ordinal_type. The cast to int must succeed.
307  fillLevel = as<int> (params.get<global_ordinal_type> ("fact: iluk level-of-fill"));
308  gotFillLevel = true;
309  }
310  catch (InvalidParameterType&) {
311  // Try the next type.
312  }
313  // Don't catch InvalidParameterName here; we've already done that above.
314  }
315 
316  if (! gotFillLevel) {
317  try {
318  // Try magnitude_type, for compatibility with ILUT.
319  // The cast from magnitude_type to int must succeed.
320  fillLevel = as<int> (params.get<magnitude_type> ("fact: iluk level-of-fill"));
321  gotFillLevel = true;
322  }
323  catch (InvalidParameterType&) {
324  // Try the next type.
325  }
326  // Don't catch InvalidParameterName here; we've already done that above.
327  }
328 
329  if (! gotFillLevel) {
330  try {
331  // Try double, for compatibility with ILUT.
332  // The cast from double to int must succeed.
333  fillLevel = as<int> (params.get<double> ("fact: iluk level-of-fill"));
334  gotFillLevel = true;
335  }
336  catch (InvalidParameterType& e) {
337  // We're out of options. The user gave us the parameter, but it
338  // doesn't have the right type. The best thing for us to do in
339  // that case is to throw, telling the user to use the right
340  // type.
341  throw e;
342  }
343  // Don't catch InvalidParameterName here; we've already done that above.
344  }
345 
347  ! gotFillLevel,
348  std::logic_error,
349  "Ifpack2::RILUK::setParameters: We should never get here! "
350  "The method should either have read the \"fact: iluk level-of-fill\" "
351  "parameter by this point, or have thrown an exception. "
352  "Please let the Ifpack2 developers know about this bug.");
353 
354  //
355  // For the other parameters, we prefer magnitude_type, but allow
356  // double for backwards compatibility.
357  //
358 
359  try {
360  absThresh = params.get<magnitude_type> ("fact: absolute threshold");
361  }
362  catch (InvalidParameterType&) {
363  // Try double, for backwards compatibility.
364  // The cast from double to magnitude_type must succeed.
365  absThresh = as<magnitude_type> (params.get<double> ("fact: absolute threshold"));
366  }
367  catch (InvalidParameterName&) {
368  // Accept the default value.
369  }
370 
371  try {
372  relThresh = params.get<magnitude_type> ("fact: relative threshold");
373  }
374  catch (InvalidParameterType&) {
375  // Try double, for backwards compatibility.
376  // The cast from double to magnitude_type must succeed.
377  relThresh = as<magnitude_type> (params.get<double> ("fact: relative threshold"));
378  }
379  catch (InvalidParameterName&) {
380  // Accept the default value.
381  }
382 
383  try {
384  relaxValue = params.get<magnitude_type> ("fact: relax value");
385  }
386  catch (InvalidParameterType&) {
387  // Try double, for backwards compatibility.
388  // The cast from double to magnitude_type must succeed.
389  relaxValue = as<magnitude_type> (params.get<double> ("fact: relax value"));
390  }
391  catch (InvalidParameterName&) {
392  // Accept the default value.
393  }
394 
395  // Forward to trisolvers.
396  L_solver_->setParameters(params);
397  U_solver_->setParameters(params);
398 
399  // "Commit" the values only after validating all of them. This
400  // ensures that there are no side effects if this routine throws an
401  // exception.
402 
403  LevelOfFill_ = fillLevel;
404 
405  // mfh 28 Nov 2012: The previous code would not assign Athresh_,
406  // Rthresh_, or RelaxValue_, if the read-in value was -1. I don't
407  // know if keeping this behavior is correct, but I'll keep it just
408  // so as not to change previous behavior.
409 
410  if (absThresh != -STM::one ()) {
411  Athresh_ = absThresh;
412  }
413  if (relThresh != -STM::one ()) {
414  Rthresh_ = relThresh;
415  }
416  if (relaxValue != -STM::one ()) {
417  RelaxValue_ = relaxValue;
418  }
419 }
420 
421 
422 template<class MatrixType>
425  return Teuchos::rcp_implicit_cast<const row_matrix_type> (A_);
426 }
427 
428 
429 template<class MatrixType>
432  return Teuchos::rcp_dynamic_cast<const crs_matrix_type> (A_, true);
433 }
434 
435 
436 template<class MatrixType>
439 {
440  using Teuchos::RCP;
441  using Teuchos::rcp;
442  using Teuchos::rcp_dynamic_cast;
443  using Teuchos::rcp_implicit_cast;
444 
445  // If A_'s communicator only has one process, or if its column and
446  // row Maps are the same, then it is already local, so use it
447  // directly.
448  if (A->getRowMap ()->getComm ()->getSize () == 1 ||
449  A->getRowMap ()->isSameAs (* (A->getColMap ()))) {
450  return A;
451  }
452 
453  // If A_ is already a LocalFilter, then use it directly. This
454  // should be the case if RILUK is being used through
455  // AdditiveSchwarz, for example.
456  RCP<const LocalFilter<row_matrix_type> > A_lf_r =
457  rcp_dynamic_cast<const LocalFilter<row_matrix_type> > (A);
458  if (! A_lf_r.is_null ()) {
459  return rcp_implicit_cast<const row_matrix_type> (A_lf_r);
460  }
461  else {
462  // A_'s communicator has more than one process, its row Map and
463  // its column Map differ, and A_ is not a LocalFilter. Thus, we
464  // have to wrap it in a LocalFilter.
465  return rcp (new LocalFilter<row_matrix_type> (A));
466  }
467 }
468 
469 
470 template<class MatrixType>
472 {
473  using Teuchos::RCP;
474  using Teuchos::rcp;
475  using Teuchos::rcp_const_cast;
476  using Teuchos::rcp_dynamic_cast;
477  using Teuchos::rcp_implicit_cast;
478  typedef Tpetra::CrsGraph<local_ordinal_type,
480  node_type> crs_graph_type;
481  const char prefix[] = "Ifpack2::RILUK::initialize: ";
482 
484  (A_.is_null (), std::runtime_error, prefix << "The matrix is null. Please "
485  "call setMatrix() with a nonnull input before calling this method.");
487  (! A_->isFillComplete (), std::runtime_error, prefix << "The matrix is not "
488  "fill complete. You may not invoke initialize() or compute() with this "
489  "matrix until the matrix is fill complete. If your matrix is a "
490  "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
491  "range Maps, if appropriate) before calling this method.");
492 
493  Teuchos::Time timer ("RILUK::initialize");
494  { // Start timing
495  Teuchos::TimeMonitor timeMon (timer);
496 
497  // Calling initialize() means that the user asserts that the graph
498  // of the sparse matrix may have changed. We must not just reuse
499  // the previous graph in that case.
500  //
501  // Regarding setting isAllocated_ to false: Eventually, we may want
502  // some kind of clever memory reuse strategy, but it's always
503  // correct just to blow everything away and start over.
504  isInitialized_ = false;
505  isAllocated_ = false;
506  isComputed_ = false;
507  Graph_ = Teuchos::null;
508 
509  A_local_ = makeLocalFilter (A_);
511  A_local_.is_null (), std::logic_error, "Ifpack2::RILUK::initialize: "
512  "makeLocalFilter returned null; it failed to compute A_local. "
513  "Please report this bug to the Ifpack2 developers.");
514 
515  // FIXME (mfh 24 Jan 2014, 26 Mar 2014) It would be more efficient
516  // to rewrite RILUK so that it works with any RowMatrix input, not
517  // just CrsMatrix. (That would require rewriting IlukGraph to
518  // handle a Tpetra::RowGraph.) However, to make it work for now,
519  // we just copy the input matrix if it's not a CrsMatrix.
520  {
521  RCP<const crs_matrix_type> A_local_crs =
522  rcp_dynamic_cast<const crs_matrix_type> (A_local_);
523  if (A_local_crs.is_null ()) {
524  // FIXME (mfh 24 Jan 2014) It would be smarter to count up the
525  // number of elements in each row of A_local, so that we can
526  // create A_local_crs_nc using static profile. The code below is
527  // correct but potentially slow.
528  RCP<crs_matrix_type> A_local_crs_nc =
529  rcp (new crs_matrix_type (A_local_->getRowMap (),
530  A_local_->getColMap (), 0));
531  // FIXME (mfh 24 Jan 2014) This Import approach will only work
532  // if A_ has a one-to-one row Map. This is generally the case
533  // with matrices given to Ifpack2.
534  //
535  // Source and destination Maps are the same in this case.
536  // That way, the Import just implements a copy.
537  typedef Tpetra::Import<local_ordinal_type, global_ordinal_type,
538  node_type> import_type;
539  import_type import (A_local_->getRowMap (), A_local_->getRowMap ());
540  A_local_crs_nc->doImport (*A_local_, import, Tpetra::REPLACE);
541  A_local_crs_nc->fillComplete (A_local_->getDomainMap (), A_local_->getRangeMap ());
542  A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
543  }
544  Graph_ = rcp (new Ifpack2::IlukGraph<crs_graph_type> (A_local_crs->getCrsGraph (),
545  LevelOfFill_, 0));
546  }
547 
548  Graph_->initialize ();
549  allocate_L_and_U ();
550  checkOrderingConsistency (*A_local_);
551  L_solver_->setMatrix (L_);
552  L_solver_->initialize ();
553  U_solver_->setMatrix (U_);
554  U_solver_->initialize ();
555  // Do not call initAllValues. compute() always calls initAllValues to
556  // fill L and U with possibly new numbers. initialize() is concerned
557  // only with the nonzero pattern.
558 
559  } // Stop timing
560 
561  isInitialized_ = true;
562  ++numInitialize_;
563  initializeTime_ += timer.totalElapsedTime ();
564 }
565 
566 template<class MatrixType>
567 void
569 checkOrderingConsistency (const row_matrix_type& A)
570 {
571  // First check that the local row map ordering is the same as the local portion of the column map.
572  // The extraction of the strictly lower/upper parts of A, as well as the factorization,
573  // implicitly assume that this is the case.
574  Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getNodeElementList();
575  Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getNodeElementList();
576  bool gidsAreConsistentlyOrdered=true;
577  global_ordinal_type indexOfInconsistentGID=0;
578  for (global_ordinal_type i=0; i<rowGIDs.size(); ++i) {
579  if (rowGIDs[i] != colGIDs[i]) {
580  gidsAreConsistentlyOrdered=false;
581  indexOfInconsistentGID=i;
582  break;
583  }
584  }
585  TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==false, std::runtime_error,
586  "The ordering of the local GIDs in the row and column maps is not the same"
587  << std::endl << "at index " << indexOfInconsistentGID
588  << ". Consistency is required, as all calculations are done with"
589  << std::endl << "local indexing.");
590 }
591 
592 template<class MatrixType>
593 void
594 RILUK<MatrixType>::
595 initAllValues (const row_matrix_type& A)
596 {
597  using Teuchos::ArrayRCP;
598  using Teuchos::Comm;
599  using Teuchos::ptr;
600  using Teuchos::RCP;
601  using Teuchos::REDUCE_SUM;
602  using Teuchos::reduceAll;
603  typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type;
604 
605  size_t NumIn = 0, NumL = 0, NumU = 0;
606  bool DiagFound = false;
607  size_t NumNonzeroDiags = 0;
608  size_t MaxNumEntries = A.getGlobalMaxNumRowEntries();
609 
610  // Allocate temporary space for extracting the strictly
611  // lower and upper parts of the matrix A.
612  Teuchos::Array<local_ordinal_type> InI(MaxNumEntries);
613  Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
614  Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
615  Teuchos::Array<scalar_type> InV(MaxNumEntries);
616  Teuchos::Array<scalar_type> LV(MaxNumEntries);
617  Teuchos::Array<scalar_type> UV(MaxNumEntries);
618 
619  // Check if values should be inserted or replaced
620  const bool ReplaceValues = L_->isStaticGraph () || L_->isLocallyIndexed ();
621 
622  L_->resumeFill ();
623  U_->resumeFill ();
624  if (ReplaceValues) {
625  L_->setAllToScalar (STS::zero ()); // Zero out L and U matrices
626  U_->setAllToScalar (STS::zero ());
627  }
628 
629  D_->putScalar (STS::zero ()); // Set diagonal values to zero
630  ArrayRCP<scalar_type> DV = D_->get1dViewNonConst (); // Get view of diagonal
631 
632  RCP<const map_type> rowMap = L_->getRowMap ();
633 
634  // First we copy the user's matrix into L and U, regardless of fill level.
635  // It is important to note that L and U are populated using local indices.
636  // This means that if the row map GIDs are not monotonically increasing
637  // (i.e., permuted or gappy), then the strictly lower (upper) part of the
638  // matrix is not the one that you would get if you based L (U) on GIDs.
639  // This is ok, as the *order* of the GIDs in the rowmap is a better
640  // expression of the user's intent than the GIDs themselves.
641 
642  Teuchos::ArrayView<const global_ordinal_type> nodeGIDs = rowMap->getNodeElementList();
643  for (size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
644  local_ordinal_type local_row = myRow;
645 
646  //TODO JJH 4April2014 An optimization is to use getLocalRowView. Not all matrices support this,
647  // we'd need to check via the Tpetra::RowMatrix method supportsRowViews().
648  A.getLocalRowCopy (local_row, InI(), InV(), NumIn); // Get Values and Indices
649 
650  // Split into L and U (we don't assume that indices are ordered).
651 
652  NumL = 0;
653  NumU = 0;
654  DiagFound = false;
655 
656  for (size_t j = 0; j < NumIn; ++j) {
657  const local_ordinal_type k = InI[j];
658 
659  if (k == local_row) {
660  DiagFound = true;
661  // Store perturbed diagonal in Tpetra::Vector D_
662  DV[local_row] += Rthresh_ * InV[j] + IFPACK2_SGN(InV[j]) * Athresh_;
663  }
664  else if (k < 0) { // Out of range
666  true, std::runtime_error, "Ifpack2::RILUK::initAllValues: current "
667  "GID k = " << k << " < 0. I'm not sure why this is an error; it is "
668  "probably an artifact of the undocumented assumptions of the "
669  "original implementation (likely copied and pasted from Ifpack). "
670  "Nevertheless, the code I found here insisted on this being an error "
671  "state, so I will throw an exception here.");
672  }
673  else if (k < local_row) {
674  LI[NumL] = k;
675  LV[NumL] = InV[j];
676  NumL++;
677  }
678  else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
679  UI[NumU] = k;
680  UV[NumU] = InV[j];
681  NumU++;
682  }
683  }
684 
685  // Check in things for this row of L and U
686 
687  if (DiagFound) {
688  ++NumNonzeroDiags;
689  } else {
690  DV[local_row] = Athresh_;
691  }
692 
693  if (NumL) {
694  if (ReplaceValues) {
695  L_->replaceLocalValues(local_row, LI(0, NumL), LV(0,NumL));
696  } else {
697  //FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
698  //FIXME in this row in the column locations corresponding to UI.
699  L_->insertLocalValues(local_row, LI(0,NumL), LV(0,NumL));
700  }
701  }
702 
703  if (NumU) {
704  if (ReplaceValues) {
705  U_->replaceLocalValues(local_row, UI(0,NumU), UV(0,NumU));
706  } else {
707  //FIXME JJH 24April2014 Is this correct? I believe this case is when there aren't already values
708  //FIXME in this row in the column locations corresponding to UI.
709  U_->insertLocalValues(local_row, UI(0,NumU), UV(0,NumU));
710  }
711  }
712  }
713 
714  // At this point L and U have the values of A in the structure of L
715  // and U, and diagonal vector D
716 
717  isInitialized_ = true;
718 }
719 
720 
721 template<class MatrixType>
723 {
724  using Teuchos::rcp;
725  const char prefix[] = "Ifpack2::RILUK::compute: ";
726 
727  // initialize() checks this too, but it's easier for users if the
728  // error shows them the name of the method that they actually
729  // called, rather than the name of some internally called method.
731  (A_.is_null (), std::runtime_error, prefix << "The matrix is null. Please "
732  "call setMatrix() with a nonnull input before calling this method.");
734  (! A_->isFillComplete (), std::runtime_error, prefix << "The matrix is not "
735  "fill complete. You may not invoke initialize() or compute() with this "
736  "matrix until the matrix is fill complete. If your matrix is a "
737  "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
738  "range Maps, if appropriate) before calling this method.");
739 
740  if (! isInitialized ()) {
741  initialize (); // Don't count this in the compute() time
742  }
743 
744  Teuchos::Time timer ("RILUK::compute");
745 
746  // Start timing
747  Teuchos::TimeMonitor timeMon (timer);
748 
749  isComputed_ = false;
750 
751  // Fill L and U with numbers. This supports nonzero pattern reuse by calling
752  // initialize() once and then compute() multiple times.
753  initAllValues (*A_local_);
754 
755  // MinMachNum should be officially defined, for now pick something a little
756  // bigger than IEEE underflow value
757 
758  const scalar_type MinDiagonalValue = STS::rmin ();
759  const scalar_type MaxDiagonalValue = STS::one () / MinDiagonalValue;
760 
761  size_t NumIn, NumL, NumU;
762 
763  // Get Maximum Row length
764  const size_t MaxNumEntries =
765  L_->getNodeMaxNumRowEntries () + U_->getNodeMaxNumRowEntries () + 1;
766 
767  Teuchos::Array<local_ordinal_type> InI(MaxNumEntries); // Allocate temp space
768  Teuchos::Array<scalar_type> InV(MaxNumEntries);
769  size_t num_cols = U_->getColMap()->getNodeNumElements();
770  Teuchos::Array<int> colflag(num_cols);
771 
772  Teuchos::ArrayRCP<scalar_type> DV = D_->get1dViewNonConst(); // Get view of diagonal
773 
774  // Now start the factorization.
775 
776  // Need some integer workspace and pointers
777  size_t NumUU;
780  for (size_t j = 0; j < num_cols; ++j) {
781  colflag[j] = -1;
782  }
783 
784  for (size_t i = 0; i < L_->getNodeNumRows (); ++i) {
785  local_ordinal_type local_row = i;
786 
787  // Fill InV, InI with current row of L, D and U combined
788 
789  NumIn = MaxNumEntries;
790  L_->getLocalRowCopy (local_row, InI (), InV (), NumL);
791 
792  InV[NumL] = DV[i]; // Put in diagonal
793  InI[NumL] = local_row;
794 
795  U_->getLocalRowCopy (local_row, InI (NumL+1, MaxNumEntries-NumL-1),
796  InV (NumL+1, MaxNumEntries-NumL-1), NumU);
797  NumIn = NumL+NumU+1;
798 
799  // Set column flags
800  for (size_t j = 0; j < NumIn; ++j) {
801  colflag[InI[j]] = j;
802  }
803 
804  scalar_type diagmod = STS::zero (); // Off-diagonal accumulator
805 
806  for (size_t jj = 0; jj < NumL; ++jj) {
807  local_ordinal_type j = InI[jj];
808  scalar_type multiplier = InV[jj]; // current_mults++;
809 
810  InV[jj] *= DV[j];
811 
812  U_->getLocalRowView(j, UUI, UUV); // View of row above
813  NumUU = UUI.size();
814 
815  if (RelaxValue_ == STM::zero ()) {
816  for (size_t k = 0; k < NumUU; ++k) {
817  const int kk = colflag[UUI[k]];
818  // FIXME (mfh 23 Dec 2013) Wait a second, we just set
819  // colflag above using size_t (which is generally unsigned),
820  // but now we're querying it using int (which is signed).
821  if (kk > -1) {
822  InV[kk] -= multiplier * UUV[k];
823  }
824  }
825  }
826  else {
827  for (size_t k = 0; k < NumUU; ++k) {
828  // FIXME (mfh 23 Dec 2013) Wait a second, we just set
829  // colflag above using size_t (which is generally unsigned),
830  // but now we're querying it using int (which is signed).
831  const int kk = colflag[UUI[k]];
832  if (kk > -1) {
833  InV[kk] -= multiplier*UUV[k];
834  }
835  else {
836  diagmod -= multiplier*UUV[k];
837  }
838  }
839  }
840  }
841  if (NumL) {
842  // Replace current row of L
843  L_->replaceLocalValues (local_row, InI (0, NumL), InV (0, NumL));
844  }
845 
846  DV[i] = InV[NumL]; // Extract Diagonal value
847 
848  if (RelaxValue_ != STM::zero ()) {
849  DV[i] += RelaxValue_*diagmod; // Add off diagonal modifications
850  }
851 
852  if (STS::magnitude (DV[i]) > STS::magnitude (MaxDiagonalValue)) {
853  if (STS::real (DV[i]) < STM::zero ()) {
854  DV[i] = -MinDiagonalValue;
855  }
856  else {
857  DV[i] = MinDiagonalValue;
858  }
859  }
860  else {
861  DV[i] = STS::one () / DV[i]; // Invert diagonal value
862  }
863 
864  for (size_t j = 0; j < NumU; ++j) {
865  InV[NumL+1+j] *= DV[i]; // Scale U by inverse of diagonal
866  }
867 
868  if (NumU) {
869  // Replace current row of L and U
870  U_->replaceLocalValues (local_row, InI (NumL+1, NumU), InV (NumL+1, NumU));
871  }
872 
873  // Reset column flags
874  for (size_t j = 0; j < NumIn; ++j) {
875  colflag[InI[j]] = -1;
876  }
877  }
878 
879  // The domain of L and the range of U are exactly their own row maps
880  // (there is no communication). The domain of U and the range of L
881  // must be the same as those of the original matrix, However if the
882  // original matrix is a VbrMatrix, these two latter maps are
883  // translation from a block map to a point map.
884  // FIXME (mfh 23 Dec 2013) Do we know that the column Map of L_ is
885  // always one-to-one?
886  L_->fillComplete (L_->getColMap (), A_local_->getRangeMap ());
887  U_->fillComplete (A_local_->getDomainMap (), U_->getRowMap ());
888 
889  // If L_solver_ or U_solver store modified factors internally, we need to reset those
890  L_solver_->setMatrix (L_);
891  L_solver_->compute ();
892  U_solver_->setMatrix (U_);
893  U_solver_->compute ();
894 
895  isComputed_ = true;
896  ++numCompute_;
897  computeTime_ += timer.totalElapsedTime ();
898 }
899 
900 
901 template<class MatrixType>
902 void
904 apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
905  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
906  Teuchos::ETransp mode,
907  scalar_type alpha,
908  scalar_type beta) const
909 {
910  using Teuchos::RCP;
911  using Teuchos::rcpFromRef;
912 
914  A_.is_null (), std::runtime_error, "Ifpack2::RILUK::apply: The matrix is "
915  "null. Please call setMatrix() with a nonnull input, then initialize() "
916  "and compute(), before calling this method.");
918  ! isComputed (), std::runtime_error,
919  "Ifpack2::RILUK::apply: If you have not yet called compute(), "
920  "you must call compute() before calling this method.");
921  TEUCHOS_TEST_FOR_EXCEPTION(
922  X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
923  "Ifpack2::RILUK::apply: X and Y do not have the same number of columns. "
924  "X.getNumVectors() = " << X.getNumVectors ()
925  << " != Y.getNumVectors() = " << Y.getNumVectors () << ".");
926  TEUCHOS_TEST_FOR_EXCEPTION(
927  STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
928  "Ifpack2::RILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
929  "complex Scalar type. Please talk to the Ifpack2 developers to get this "
930  "fixed. There is a FIXME in this file about this very issue.");
931 #ifdef HAVE_IFPACK2_DEBUG
932  {
933  const magnitude_type D_nrm1 = D_->norm1 ();
934  TEUCHOS_TEST_FOR_EXCEPTION( STM::isnaninf (D_nrm1), std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the stored diagonal is NaN or Inf.");
935  Teuchos::Array<magnitude_type> norms (X.getNumVectors ());
936  X.norm1 (norms ());
937  bool good = true;
938  for (size_t j = 0; j < X.getNumVectors (); ++j) {
939  if (STM::isnaninf (norms[j])) {
940  good = false;
941  break;
942  }
943  }
944  TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the input X is NaN or Inf.");
945  }
946 #endif // HAVE_IFPACK2_DEBUG
947 
948  const scalar_type one = STS::one ();
949  const scalar_type zero = STS::zero ();
950 
951  Teuchos::Time timer ("RILUK::apply");
952  { // Start timing
953  Teuchos::TimeMonitor timeMon (timer);
954  if (alpha == one && beta == zero) {
955  if (mode == Teuchos::NO_TRANS) { // Solve L (D (U Y)) = X for Y.
956  // Start by solving L Y = X for Y.
957  L_solver_->apply (X, Y, mode);
958 
959  // Solve D Y = Y. The operation lets us do this in place in Y, so we can
960  // write "solve D Y = Y for Y."
961  Y.elementWiseMultiply (one, *D_, Y, zero);
962 
963  U_solver_->apply (Y, Y, mode); // Solve U Y = Y.
964  }
965  else { // Solve U^P (D^P (L^P Y)) = X for Y (where P is * or T).
966 
967  // Start by solving U^P Y = X for Y.
968  U_solver_->apply (X, Y, mode);
969 
970  // Solve D^P Y = Y.
971  //
972  // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we
973  // need to do an elementwise multiply with the conjugate of
974  // D_, not just with D_ itself.
975  Y.elementWiseMultiply (one, *D_, Y, zero);
976 
977  L_solver_->apply (Y, Y, mode); // Solve L^P Y = Y.
978  }
979  }
980  else { // alpha != 1 or beta != 0
981  if (alpha == zero) {
982  // The special case for beta == 0 ensures that if Y contains Inf
983  // or NaN values, we replace them with 0 (following BLAS
984  // convention), rather than multiplying them by 0 to get NaN.
985  if (beta == zero) {
986  Y.putScalar (zero);
987  } else {
988  Y.scale (beta);
989  }
990  } else { // alpha != zero
991  MV Y_tmp (Y.getMap (), Y.getNumVectors ());
992  apply (X, Y_tmp, mode);
993  Y.update (alpha, Y_tmp, beta);
994  }
995  }
996  }//end timing
997 
998 #ifdef HAVE_IFPACK2_DEBUG
999  {
1000  Teuchos::Array<magnitude_type> norms (Y.getNumVectors ());
1001  Y.norm1 (norms ());
1002  bool good = true;
1003  for (size_t j = 0; j < Y.getNumVectors (); ++j) {
1004  if (STM::isnaninf (norms[j])) {
1005  good = false;
1006  break;
1007  }
1008  }
1009  TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error, "Ifpack2::RILUK::apply: The 1-norm of the output Y is NaN or Inf.");
1010  }
1011 #endif // HAVE_IFPACK2_DEBUG
1012 
1013  ++numApply_;
1014  applyTime_ = timer.totalElapsedTime ();
1015 }
1016 
1017 
1018 template<class MatrixType>
1020 multiply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
1021  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
1022  const Teuchos::ETransp mode) const
1023 {
1024  const scalar_type zero = STS::zero ();
1025  const scalar_type one = STS::one ();
1026 
1027  if (mode != Teuchos::NO_TRANS) {
1028  U_->apply (X, Y, mode); //
1029  Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1030 
1031  // FIXME (mfh 24 Jan 2014) If mode = Teuchos::CONJ_TRANS, we need
1032  // to do an elementwise multiply with the conjugate of D_, not
1033  // just with D_ itself.
1034  Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1035 
1036  MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y
1037  L_->apply (Y_tmp, Y, mode);
1038  Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1039  }
1040  else {
1041  L_->apply (X, Y, mode);
1042  Y.update (one, X, one); // Y = Y + X (account for implicit unit diagonal)
1043  Y.elementWiseMultiply (one, *D_, Y, zero); // y = D*y (D_ has inverse of diagonal)
1044  MV Y_tmp (Y, Teuchos::Copy); // Need a temp copy of Y1
1045  U_->apply (Y_tmp, Y, mode);
1046  Y.update (one, Y_tmp, one); // (account for implicit unit diagonal)
1047  }
1048 }
1049 
1050 
1051 template<class MatrixType>
1053 {
1054  std::ostringstream os;
1055 
1056  // Output is a valid YAML dictionary in flow style. If you don't
1057  // like everything on a single line, you should call describe()
1058  // instead.
1059  os << "\"Ifpack2::RILUK\": {";
1060  os << "Initialized: " << (isInitialized () ? "true" : "false") << ", "
1061  << "Computed: " << (isComputed () ? "true" : "false") << ", ";
1062 
1063  os << "Level-of-fill: " << getLevelOfFill() << ", ";
1064 
1065  if (A_.is_null ()) {
1066  os << "Matrix: null";
1067  }
1068  else {
1069  os << "Global matrix dimensions: ["
1070  << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"
1071  << ", Global nnz: " << A_->getGlobalNumEntries();
1072  }
1073 
1074  if (! L_solver_.is_null ()) os << ", " << L_solver_->description ();
1075  if (! U_solver_.is_null ()) os << ", " << U_solver_->description ();
1076 
1077  os << "}";
1078  return os.str ();
1079 }
1080 
1081 } // namespace Ifpack2
1082 
1083 #define IFPACK2_RILUK_INSTANT(S,LO,GO,N) \
1084  template class Ifpack2::RILUK< Tpetra::RowMatrix<S, LO, GO, N> >;
1085 
1086 #endif
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_RILUK_def.hpp:223
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:261
T & get(const std::string &name, T def_value)
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > & getD() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:159
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:471
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_RILUK_decl.hpp:273
size_type size() const
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_RILUK_def.hpp:106
std::string description() const
A one-line description of this object.
Definition: Ifpack2_RILUK_def.hpp:1052
void setParameters(const Teuchos::ParameterList &params)
Definition: Ifpack2_RILUK_def.hpp:270
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:243
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_RILUK_def.hpp:186
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition: Ifpack2_RILUK_def.hpp:431
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition: Ifpack2_RILUK_decl.hpp:276
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:173
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:264
T * getRawPtr() const
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_RILUK_def.hpp:204
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:258
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_RILUK_def.hpp:95
&quot;Preconditioner&quot; that solves local sparse triangular systems.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:77
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:722
TEUCHOS_DEPRECATED void reduceAll(const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, Packet *globalReduct)
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 (inverse of the) incomplete factorization to X, resulting in Y.
Definition: Ifpack2_RILUK_def.hpp:904
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_RILUK_decl.hpp:267
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition: Ifpack2_RILUK_def.hpp:424
TypeTo as(const TypeFrom &t)
double totalElapsedTime(bool readCurrentTime=false) const
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:142
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:255