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