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