Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_ILUT_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_ILUT_DEF_HPP
11#define IFPACK2_ILUT_DEF_HPP
12
13#include <type_traits>
14#include "Kokkos_StaticCrsGraph.hpp"
15#include "Teuchos_TypeNameTraits.hpp"
16#include "Teuchos_StandardParameterEntryValidators.hpp"
17#include "Teuchos_Time.hpp"
18#include "Tpetra_CrsMatrix.hpp"
19#include "KokkosSparse_par_ilut.hpp"
20
21#include "Ifpack2_Heap.hpp"
22#include "Ifpack2_LocalFilter.hpp"
23#include "Ifpack2_LocalSparseTriangularSolver.hpp"
24#include "Ifpack2_Parameters.hpp"
25#include "Ifpack2_Details_getParamTryingTypes.hpp"
26
27namespace Ifpack2 {
28
29 namespace {
30
31 struct IlutImplType {
32 enum Enum {
33 Serial,
34 PAR_ILUT
35 };
36
37 static void loadPLTypeOption (Teuchos::Array<std::string>& type_strs, Teuchos::Array<Enum>& type_enums) {
38 type_strs.resize(2);
39 type_strs[0] = "serial";
40 type_strs[1] = "par_ilut";
41 type_enums.resize(2);
42 type_enums[0] = Serial;
43 type_enums[1] = PAR_ILUT;
44 }
45 };
46
47
72 template<class ScalarType>
73 inline typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
74 ilutDefaultDropTolerance () {
75 typedef Teuchos::ScalarTraits<ScalarType> STS;
76 typedef typename STS::magnitudeType magnitude_type;
77 typedef Teuchos::ScalarTraits<magnitude_type> STM;
78
79 // 1/2. Hopefully this can be represented in magnitude_type.
80 const magnitude_type oneHalf = STM::one() / (STM::one() + STM::one());
81
82 // The min ensures that in case magnitude_type has very low
83 // precision, we'll at least get some value strictly less than
84 // one.
85 return std::min (static_cast<magnitude_type> (1000) * STS::magnitude (STS::eps ()), oneHalf);
86 }
87
88 // Full specialization for ScalarType = double.
89 // This specialization preserves ILUT's previous default behavior.
90 template<>
91 inline Teuchos::ScalarTraits<double>::magnitudeType
92 ilutDefaultDropTolerance<double> () {
93 return 1e-12;
94 }
95
96 } // namespace (anonymous)
97
98
99template <class MatrixType>
100ILUT<MatrixType>::ILUT (const Teuchos::RCP<const row_matrix_type>& A) :
101 A_ (A),
102 Athresh_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
103 Rthresh_ (Teuchos::ScalarTraits<magnitude_type>::one ()),
104 RelaxValue_ (Teuchos::ScalarTraits<magnitude_type>::zero ()),
105 LevelOfFill_ (1.0),
106 DropTolerance_ (ilutDefaultDropTolerance<scalar_type> ()),
107 par_ilut_options_{1, 0., -1, -1, 0.75, false},
108 InitializeTime_ (0.0),
109 ComputeTime_ (0.0),
110 ApplyTime_ (0.0),
111 NumInitialize_ (0),
112 NumCompute_ (0),
113 NumApply_ (0),
114 IsInitialized_ (false),
115 IsComputed_ (false),
116 useKokkosKernelsParILUT_(false)
117
118{
119 allocateSolvers();
120}
121
122template<class MatrixType>
123void ILUT<MatrixType>::allocateSolvers ()
124{
125 L_solver_ = Teuchos::rcp (new LocalSparseTriangularSolver<row_matrix_type> ());
126 L_solver_->setObjectLabel("lower");
127 U_solver_ = Teuchos::rcp (new LocalSparseTriangularSolver<row_matrix_type> ());
128 U_solver_->setObjectLabel("upper");
129}
130
131template <class MatrixType>
132void ILUT<MatrixType>::setParameters (const Teuchos::ParameterList& params)
133{
134 using Ifpack2::Details::getParamTryingTypes;
135 const char prefix[] = "Ifpack2::ILUT: ";
136
137 // Don't actually change the instance variables until we've checked
138 // all parameters. This ensures that setParameters satisfies the
139 // strong exception guarantee (i.e., is transactional).
140
141 // Parsing implementation type
142 IlutImplType::Enum ilutimplType = IlutImplType::Serial;
143 do {
144 static const char typeName[] = "fact: type";
145
146 if ( ! params.isType<std::string>(typeName)) break;
147
148 // Map std::string <-> IlutImplType::Enum.
149 Teuchos::Array<std::string> ilutimplTypeStrs;
150 Teuchos::Array<IlutImplType::Enum> ilutimplTypeEnums;
151 IlutImplType::loadPLTypeOption (ilutimplTypeStrs, ilutimplTypeEnums);
152 Teuchos::StringToIntegralParameterEntryValidator<IlutImplType::Enum>
153 s2i(ilutimplTypeStrs (), ilutimplTypeEnums (), typeName, false);
154
155 ilutimplType = s2i.getIntegralValue(params.get<std::string>(typeName));
156 } while (0);
157
158 if (ilutimplType == IlutImplType::PAR_ILUT) {
159 this->useKokkosKernelsParILUT_ = true;
160 }
161 else {
162 this->useKokkosKernelsParILUT_ = false;
163 }
164
165 // Fill level in ILUT is a double, not a magnitude_type, because it
166 // depends on LO and GO, not on Scalar. Also, you can't cast
167 // arbitrary magnitude_type (e.g., Sacado::MP::Vector) to double.
168 double fillLevel = LevelOfFill_;
169 {
170 const std::string paramName ("fact: ilut level-of-fill");
171 TEUCHOS_TEST_FOR_EXCEPTION(
172 (params.isParameter(paramName) && this->useKokkosKernelsParILUT_), std::runtime_error,
173 "Ifpack2::ILUT: Parameter " << paramName << " is meaningless for algorithm par_ilut.");
174 getParamTryingTypes<double, double, float>
175 (fillLevel, params, paramName, prefix);
176 TEUCHOS_TEST_FOR_EXCEPTION
177 (fillLevel < 1.0, std::runtime_error,
178 "Ifpack2::ILUT: The \"" << paramName << "\" parameter must be >= "
179 "1.0, but you set it to " << fillLevel << ". For ILUT, the fill level "
180 "means something different than it does for ILU(k). ILU(0) produces "
181 "factors with the same sparsity structure as the input matrix A. For "
182 "ILUT, level-of-fill = 1.0 will produce factors with nonzeros matching "
183 "the sparsity structure of A. level-of-fill > 1.0 allows for additional "
184 "fill-in.");
185 }
186
187 magnitude_type absThresh = Athresh_;
188 {
189 const std::string paramName ("fact: absolute threshold");
190 getParamTryingTypes<magnitude_type, magnitude_type, double>
191 (absThresh, params, paramName, prefix);
192 }
193
194 magnitude_type relThresh = Rthresh_;
195 {
196 const std::string paramName ("fact: relative threshold");
197 getParamTryingTypes<magnitude_type, magnitude_type, double>
198 (relThresh, params, paramName, prefix);
199 }
200
201 magnitude_type relaxValue = RelaxValue_;
202 {
203 const std::string paramName ("fact: relax value");
204 getParamTryingTypes<magnitude_type, magnitude_type, double>
205 (relaxValue, params, paramName, prefix);
206 }
207
208 magnitude_type dropTol = DropTolerance_;
209 {
210 const std::string paramName ("fact: drop tolerance");
211 getParamTryingTypes<magnitude_type, magnitude_type, double>
212 (dropTol, params, paramName, prefix);
213 }
214
215 int par_ilut_max_iter=20;
216 magnitude_type par_ilut_residual_norm_delta_stop=1e-2;
217 int par_ilut_team_size=0;
218 int par_ilut_vector_size=0;
219 float par_ilut_fill_in_limit=0.75;
220 bool par_ilut_verbose=false;
221 if (this->useKokkosKernelsParILUT_) {
222 par_ilut_max_iter = par_ilut_options_.max_iter;
223 par_ilut_residual_norm_delta_stop = par_ilut_options_.residual_norm_delta_stop;
224 par_ilut_team_size = par_ilut_options_.team_size;
225 par_ilut_vector_size = par_ilut_options_.vector_size;
226 par_ilut_fill_in_limit = par_ilut_options_.fill_in_limit;
227 par_ilut_verbose = par_ilut_options_.verbose;
228
229 std::string par_ilut_plist_name("parallel ILUT options");
230 if (params.isSublist(par_ilut_plist_name)) {
231 Teuchos::ParameterList const &par_ilut_plist = params.sublist(par_ilut_plist_name);
232
233 std::string paramName("maximum iterations");
234 getParamTryingTypes<int, int>(par_ilut_max_iter, par_ilut_plist, paramName, prefix);
235
236 paramName = "residual norm delta stop";
237 getParamTryingTypes<magnitude_type, magnitude_type, double>(par_ilut_residual_norm_delta_stop, par_ilut_plist, paramName, prefix);
238
239 paramName = "team size";
240 getParamTryingTypes<int, int>(par_ilut_team_size, par_ilut_plist, paramName, prefix);
241
242 paramName = "vector size";
243 getParamTryingTypes<int, int>(par_ilut_vector_size, par_ilut_plist, paramName, prefix);
244
245 paramName = "fill in limit";
246 getParamTryingTypes<float, float, double>(par_ilut_fill_in_limit, par_ilut_plist, paramName, prefix);
247
248 paramName = "verbose";
249 getParamTryingTypes<bool, bool>(par_ilut_verbose, par_ilut_plist, paramName, prefix);
250
251 } // if (params.isSublist(par_ilut_plist_name))
252
253 par_ilut_options_.max_iter = par_ilut_max_iter;
254 par_ilut_options_.residual_norm_delta_stop = par_ilut_residual_norm_delta_stop;
255 par_ilut_options_.team_size = par_ilut_team_size;
256 par_ilut_options_.vector_size = par_ilut_vector_size;
257 par_ilut_options_.fill_in_limit = par_ilut_fill_in_limit;
258 par_ilut_options_.verbose = par_ilut_verbose;
259
260 } //if (this->useKokkosKernelsParILUT_)
261
262 // Forward to trisolvers.
263 L_solver_->setParameters(params);
264 U_solver_->setParameters(params);
265
266 LevelOfFill_ = fillLevel;
267 Athresh_ = absThresh;
268 Rthresh_ = relThresh;
269 RelaxValue_ = relaxValue;
270 DropTolerance_ = dropTol;
271}
272
273
274template <class MatrixType>
275Teuchos::RCP<const Teuchos::Comm<int> >
277 TEUCHOS_TEST_FOR_EXCEPTION(
278 A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getComm: "
279 "The matrix is null. Please call setMatrix() with a nonnull input "
280 "before calling this method.");
281 return A_->getComm ();
282}
283
284
285template <class MatrixType>
286Teuchos::RCP<const typename ILUT<MatrixType>::row_matrix_type>
288 return A_;
289}
290
291
292template <class MatrixType>
293Teuchos::RCP<const typename ILUT<MatrixType>::map_type>
295{
296 TEUCHOS_TEST_FOR_EXCEPTION(
297 A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getDomainMap: "
298 "The matrix is null. Please call setMatrix() with a nonnull input "
299 "before calling this method.");
300 return A_->getDomainMap ();
301}
302
303
304template <class MatrixType>
305Teuchos::RCP<const typename ILUT<MatrixType>::map_type>
307{
308 TEUCHOS_TEST_FOR_EXCEPTION(
309 A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getRangeMap: "
310 "The matrix is null. Please call setMatrix() with a nonnull input "
311 "before calling this method.");
312 return A_->getRangeMap ();
313}
314
315
316template <class MatrixType>
318 return true;
319}
320
321
322template <class MatrixType>
324 return NumInitialize_;
325}
326
327
328template <class MatrixType>
330 return NumCompute_;
331}
332
333
334template <class MatrixType>
336 return NumApply_;
337}
338
339
340template <class MatrixType>
342 return InitializeTime_;
343}
344
345
346template<class MatrixType>
348 return ComputeTime_;
349}
350
351
352template<class MatrixType>
354 return ApplyTime_;
355}
356
357
358template<class MatrixType>
360 TEUCHOS_TEST_FOR_EXCEPTION(
361 A_.is_null (), std::runtime_error, "Ifpack2::ILUT::getNodeSmootherComplexity: "
362 "The input matrix A is null. Please call setMatrix() with a nonnull "
363 "input matrix, then call compute(), before calling this method.");
364 // ILUT methods cost roughly one apply + the nnz in the upper+lower triangles
365 return A_->getLocalNumEntries() + getLocalNumEntries();
366}
367
368
369template<class MatrixType>
371 return L_->getGlobalNumEntries () + U_->getGlobalNumEntries ();
372}
373
374
375template<class MatrixType>
377 return L_->getLocalNumEntries () + U_->getLocalNumEntries ();
378}
379
380
381template<class MatrixType>
382void ILUT<MatrixType>::setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
383{
384 if (A.getRawPtr () != A_.getRawPtr ()) {
385 // Check in serial or one-process mode if the matrix is square.
386 TEUCHOS_TEST_FOR_EXCEPTION(
387 ! A.is_null () && A->getComm ()->getSize () == 1 &&
388 A->getLocalNumRows () != A->getLocalNumCols (),
389 std::runtime_error, "Ifpack2::ILUT::setMatrix: If A's communicator only "
390 "contains one process, then A must be square. Instead, you provided a "
391 "matrix A with " << A->getLocalNumRows () << " rows and "
392 << A->getLocalNumCols () << " columns.");
393
394 // It's legal for A to be null; in that case, you may not call
395 // initialize() until calling setMatrix() with a nonnull input.
396 // Regardless, setting the matrix invalidates any previous
397 // factorization.
398 IsInitialized_ = false;
399 IsComputed_ = false;
400 A_local_ = Teuchos::null;
401
402 // The sparse triangular solvers get a triangular factor as their
403 // input matrix. The triangular factors L_ and U_ are getting
404 // reset, so we reset the solvers' matrices to null. Do that
405 // before setting L_ and U_ to null, so that latter step actually
406 // frees the factors.
407 if (! L_solver_.is_null ()) {
408 L_solver_->setMatrix (Teuchos::null);
409 }
410 if (! U_solver_.is_null ()) {
411 U_solver_->setMatrix (Teuchos::null);
412 }
413
414 L_ = Teuchos::null;
415 U_ = Teuchos::null;
416 A_ = A;
417 }
418}
419
420template <class MatrixType>
421Teuchos::RCP<const typename ILUT<MatrixType>::row_matrix_type>
422ILUT<MatrixType>::makeLocalFilter (const Teuchos::RCP<const row_matrix_type>& A)
423{
424 using Teuchos::RCP;
425 using Teuchos::rcp;
426 using Teuchos::rcp_dynamic_cast;
427 using Teuchos::rcp_implicit_cast;
428
429 // If A_'s communicator only has one process, or if its column and
430 // row Maps are the same, then it is already local, so use it
431 // directly.
432 if (A->getRowMap ()->getComm ()->getSize () == 1 ||
433 A->getRowMap ()->isSameAs (* (A->getColMap ()))) {
434 return A;
435 }
436
437 // If A_ is already a LocalFilter, then use it directly. This
438 // should be the case if RILUT is being used through
439 // AdditiveSchwarz, for example.
440 RCP<const LocalFilter<row_matrix_type> > A_lf_r =
441 rcp_dynamic_cast<const LocalFilter<row_matrix_type> > (A);
442 if (! A_lf_r.is_null ()) {
443 return rcp_implicit_cast<const row_matrix_type> (A_lf_r);
444 }
445 else {
446 // A_'s communicator has more than one process, its row Map and
447 // its column Map differ, and A_ is not a LocalFilter. Thus, we
448 // have to wrap it in a LocalFilter.
449 return rcp (new LocalFilter<row_matrix_type> (A));
450 }
451}
452
453
454template<class MatrixType>
456{
457 using Teuchos::RCP;
458 using Teuchos::Array;
459 using Teuchos::rcp_const_cast;
460 Teuchos::Time timer ("ILUT::initialize");
461 double startTime = timer.wallTime();
462 {
463 Teuchos::TimeMonitor timeMon (timer);
464
465 // Check that the matrix is nonnull.
466 TEUCHOS_TEST_FOR_EXCEPTION(
467 A_.is_null (), std::runtime_error, "Ifpack2::ILUT::initialize: "
468 "The matrix to precondition is null. Please call setMatrix() with a "
469 "nonnull input before calling this method.");
470
471 // Clear any previous computations.
472 IsInitialized_ = false;
473 IsComputed_ = false;
474 A_local_ = Teuchos::null;
475 L_ = Teuchos::null;
476 U_ = Teuchos::null;
477
478 A_local_ = makeLocalFilter(A_); // Compute the local filter.
479 TEUCHOS_TEST_FOR_EXCEPTION(
480 A_local_.is_null(), std::logic_error, "Ifpack2::RILUT::initialize: "
481 "makeLocalFilter returned null; it failed to compute A_local. "
482 "Please report this bug to the Ifpack2 developers.");
483
484 if (this->useKokkosKernelsParILUT_) {
485 this->KernelHandle_ = Teuchos::rcp(new kk_handle_type());
486 KernelHandle_->create_par_ilut_handle();
487 auto par_ilut_handle = KernelHandle_->get_par_ilut_handle();
488 par_ilut_handle->set_residual_norm_delta_stop(par_ilut_options_.residual_norm_delta_stop);
489 par_ilut_handle->set_team_size(par_ilut_options_.team_size);
490 par_ilut_handle->set_vector_size(par_ilut_options_.vector_size);
491 par_ilut_handle->set_fill_in_limit(par_ilut_options_.fill_in_limit);
492 par_ilut_handle->set_verbose(par_ilut_options_.verbose);
493 par_ilut_handle->set_async_update(false);
494
495 RCP<const crs_matrix_type> A_local_crs = Teuchos::rcp_dynamic_cast<const crs_matrix_type>(A_local_);
496 if (A_local_crs.is_null()) {
497 // the result is a host-based matrix, which is the same as what happens in RILUK
498 local_ordinal_type numRows = A_local_->getLocalNumRows();
499 Array<size_t> entriesPerRow(numRows);
500 for(local_ordinal_type i = 0; i < numRows; i++) {
501 entriesPerRow[i] = A_local_->getNumEntriesInLocalRow(i);
502 }
503 RCP<crs_matrix_type> A_local_crs_nc =
504 rcp (new crs_matrix_type (A_local_->getRowMap (),
505 A_local_->getColMap (),
506 entriesPerRow()));
507 // copy entries into A_local_crs
508 nonconst_local_inds_host_view_type indices("indices",A_local_->getLocalMaxNumRowEntries());
509 nonconst_values_host_view_type values("values",A_local_->getLocalMaxNumRowEntries());
510 for(local_ordinal_type i = 0; i < numRows; i++) {
511 size_t numEntries = 0;
512 A_local_->getLocalRowCopy(i, indices, values, numEntries);
513 A_local_crs_nc->insertLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()), indices.data());
514 }
515 A_local_crs_nc->fillComplete (A_local_->getDomainMap (), A_local_->getRangeMap ());
516 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
517 }
518 auto A_local_crs_device = A_local_crs->getLocalMatrixDevice();
519
520 //KokkosKernels requires unsigned
521 typedef typename Kokkos::View<usize_type*, array_layout, device_type> ulno_row_view_t;
522 const int NumMyRows = A_local_crs->getRowMap()->getLocalNumElements();
523 L_rowmap_ = ulno_row_view_t("L_row_map", NumMyRows + 1);
524 U_rowmap_ = ulno_row_view_t("U_row_map", NumMyRows + 1);
525 L_rowmap_orig_ = ulno_row_view_t("L_row_map_orig", NumMyRows + 1);
526 U_rowmap_orig_ = ulno_row_view_t("U_row_map_orig", NumMyRows + 1);
527
528 KokkosSparse::Experimental::par_ilut_symbolic(KernelHandle_.getRawPtr(),
529 A_local_crs_device.graph.row_map, A_local_crs_device.graph.entries,
530 L_rowmap_,
531 U_rowmap_);
532
533 Kokkos::deep_copy(L_rowmap_orig_, L_rowmap_);
534 Kokkos::deep_copy(U_rowmap_orig_, U_rowmap_);
535 }
536
537 IsInitialized_ = true;
538 ++NumInitialize_;
539 } //timer scope
540 InitializeTime_ += (timer.wallTime() - startTime);
541}
542
543
544template<typename ScalarType>
545typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
546scalar_mag (const ScalarType& s)
547{
548 return Teuchos::ScalarTraits<ScalarType>::magnitude(s);
549}
550
551
552template<class MatrixType>
554{
555 using Teuchos::Array;
556 using Teuchos::ArrayRCP;
557 using Teuchos::ArrayView;
558 using Teuchos::as;
559 using Teuchos::rcp;
560 using Teuchos::reduceAll;
561 using Teuchos::RCP;
562 using Teuchos::rcp_const_cast;
563
564 // Don't count initialization in the compute() time.
565 if (! isInitialized ()) {
566 initialize ();
567 }
568
569 Teuchos::Time timer ("ILUT::compute");
570 double startTime = timer.wallTime();
571 { // Timer scope for timing compute()
572 Teuchos::TimeMonitor timeMon (timer, true);
573
574 if (!this->useKokkosKernelsParILUT_)
575 {
576 //--------------------------------------------------------------------------
577 // Ifpack2::ILUT's serial version is a translation of the Aztec ILUT
578 // implementation. The Aztec ILUT implementation was written by Ray Tuminaro.
579 //
580 // This isn't an exact translation of the Aztec ILUT algorithm, for the
581 // following reasons:
582 // 1. Minor differences result from the fact that Aztec factors a MSR format
583 // matrix in place, while the code below factors an input CrsMatrix which
584 // remains untouched and stores the resulting factors in separate L and U
585 // CrsMatrix objects.
586 // Also, the Aztec code begins by shifting the matrix pointers back
587 // by one, and the pointer contents back by one, and then using 1-based
588 // Fortran-style indexing in the algorithm. This Ifpack2 code uses C-style
589 // 0-based indexing throughout.
590 // 2. Aztec stores the inverse of the diagonal of U. This Ifpack2 code
591 // stores the non-inverted diagonal in U.
592 // The triangular solves (in Ifpack2::ILUT::apply()) are performed by
593 // calling the Tpetra::CrsMatrix::solve method on the L and U objects, and
594 // this requires U to contain the non-inverted diagonal.
595 //
596 // ABW.
597 //--------------------------------------------------------------------------
598
599 const scalar_type zero = STS::zero ();
600 const scalar_type one = STS::one ();
601
602 const local_ordinal_type myNumRows = A_local_->getLocalNumRows ();
603
604 // If this macro is defined, files containing the L and U factors
605 // will be written. DON'T CHECK IN THE CODE WITH THIS MACRO ENABLED!!!
606 // #define IFPACK2_WRITE_ILUT_FACTORS
607#ifdef IFPACK2_WRITE_ILUT_FACTORS
608 std::ofstream ofsL("L.ifpack2_ilut.mtx", std::ios::out);
609 std::ofstream ofsU("U.ifpack2_ilut.mtx", std::ios::out);
610#endif
611
612 // Calculate how much fill will be allowed in addition to the
613 // space that corresponds to the input matrix entries.
614 double local_nnz = static_cast<double> (A_local_->getLocalNumEntries ());
615 double fill = ((getLevelOfFill () - 1.0) * local_nnz) / (2 * myNumRows);
616
617 // std::ceil gives the smallest integer larger than the argument.
618 // this may give a slightly different result than Aztec's fill value in
619 // some cases.
620 double fill_ceil=std::ceil(fill);
621
622 // Similarly to Aztec, we will allow the same amount of fill for each
623 // row, half in L and half in U.
624 size_type fillL = static_cast<size_type>(fill_ceil);
625 size_type fillU = static_cast<size_type>(fill_ceil);
626
627 Array<scalar_type> InvDiagU (myNumRows, zero);
628
629 Array<Array<local_ordinal_type> > L_tmp_idx(myNumRows);
630 Array<Array<scalar_type> > L_tmpv(myNumRows);
631 Array<Array<local_ordinal_type> > U_tmp_idx(myNumRows);
632 Array<Array<scalar_type> > U_tmpv(myNumRows);
633
634 enum { UNUSED, ORIG, FILL };
635 local_ordinal_type max_col = myNumRows;
636
637 Array<int> pattern(max_col, UNUSED);
638 Array<scalar_type> cur_row(max_col, zero);
639 Array<magnitude_type> unorm(max_col);
640 magnitude_type rownorm;
641 Array<local_ordinal_type> L_cols_heap;
642 Array<local_ordinal_type> U_cols;
643 Array<local_ordinal_type> L_vals_heap;
644 Array<local_ordinal_type> U_vals_heap;
645
646 // A comparison object which will be used to create 'heaps' of indices
647 // that are ordered according to the corresponding values in the
648 // 'cur_row' array.
649 greater_indirect<scalar_type,local_ordinal_type> vals_comp(cur_row);
650
651 // =================== //
652 // start factorization //
653 // =================== //
654 nonconst_local_inds_host_view_type ColIndicesARCP;
655 nonconst_values_host_view_type ColValuesARCP;
656 if (! A_local_->supportsRowViews ()) {
657 const size_t maxnz = A_local_->getLocalMaxNumRowEntries ();
658 Kokkos::resize(ColIndicesARCP,maxnz);
659 Kokkos::resize(ColValuesARCP,maxnz);
660 }
661
662 for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
663 local_inds_host_view_type ColIndicesA;
664 values_host_view_type ColValuesA;
665 size_t RowNnz;
666
667 if (A_local_->supportsRowViews ()) {
668 A_local_->getLocalRowView (row_i, ColIndicesA, ColValuesA);
669 RowNnz = ColIndicesA.size ();
670 }
671 else {
672 A_local_->getLocalRowCopy (row_i, ColIndicesARCP, ColValuesARCP, RowNnz);
673 ColIndicesA = Kokkos::subview(ColIndicesARCP,std::make_pair((size_t)0, RowNnz));
674 ColValuesA = Kokkos::subview(ColValuesARCP,std::make_pair((size_t)0, RowNnz));
675 }
676
677 // Always include the diagonal in the U factor. The value should get
678 // set in the next loop below.
679 U_cols.push_back(row_i);
680 cur_row[row_i] = zero;
681 pattern[row_i] = ORIG;
682
683 size_type L_cols_heaplen = 0;
684 rownorm = STM::zero ();
685 for (size_t i = 0; i < RowNnz; ++i) {
686 if (ColIndicesA[i] < myNumRows) {
687 if (ColIndicesA[i] < row_i) {
688 add_to_heap(ColIndicesA[i], L_cols_heap, L_cols_heaplen);
689 }
690 else if (ColIndicesA[i] > row_i) {
691 U_cols.push_back(ColIndicesA[i]);
692 }
693
694 cur_row[ColIndicesA[i]] = ColValuesA[i];
695 pattern[ColIndicesA[i]] = ORIG;
696 rownorm += scalar_mag(ColValuesA[i]);
697 }
698 }
699
700 // Alter the diagonal according to the absolute-threshold and
701 // relative-threshold values. If not set, those values default
702 // to zero and one respectively.
703 const magnitude_type rthresh = getRelativeThreshold();
704 const scalar_type& v = cur_row[row_i];
705 cur_row[row_i] = as<scalar_type> (getAbsoluteThreshold() * IFPACK2_SGN(v)) + rthresh*v;
706
707 size_type orig_U_len = U_cols.size();
708 RowNnz = L_cols_heap.size() + orig_U_len;
709 rownorm = getDropTolerance() * rownorm/RowNnz;
710
711 // The following while loop corresponds to the 'L30' goto's in Aztec.
712 size_type L_vals_heaplen = 0;
713 while (L_cols_heaplen > 0) {
714 local_ordinal_type row_k = L_cols_heap.front();
715
716 scalar_type multiplier = cur_row[row_k] * InvDiagU[row_k];
717 cur_row[row_k] = multiplier;
718 magnitude_type mag_mult = scalar_mag(multiplier);
719 if (mag_mult*unorm[row_k] < rownorm) {
720 pattern[row_k] = UNUSED;
721 rm_heap_root(L_cols_heap, L_cols_heaplen);
722 continue;
723 }
724 if (pattern[row_k] != ORIG) {
725 if (L_vals_heaplen < fillL) {
726 add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
727 }
728 else if (L_vals_heaplen==0 ||
729 mag_mult < scalar_mag(cur_row[L_vals_heap.front()])) {
730 pattern[row_k] = UNUSED;
731 rm_heap_root(L_cols_heap, L_cols_heaplen);
732 continue;
733 }
734 else {
735 pattern[L_vals_heap.front()] = UNUSED;
736 rm_heap_root(L_vals_heap, L_vals_heaplen, vals_comp);
737 add_to_heap(row_k, L_vals_heap, L_vals_heaplen, vals_comp);
738 }
739 }
740
741 /* Reduce current row */
742
743 ArrayView<local_ordinal_type> ColIndicesU = U_tmp_idx[row_k]();
744 ArrayView<scalar_type> ColValuesU = U_tmpv[row_k]();
745 size_type ColNnzU = ColIndicesU.size();
746
747 for(size_type j=0; j<ColNnzU; ++j) {
748 if (ColIndicesU[j] > row_k) {
749 scalar_type tmp = multiplier * ColValuesU[j];
750 local_ordinal_type col_j = ColIndicesU[j];
751 if (pattern[col_j] != UNUSED) {
752 cur_row[col_j] -= tmp;
753 }
754 else if (scalar_mag(tmp) > rownorm) {
755 cur_row[col_j] = -tmp;
756 pattern[col_j] = FILL;
757 if (col_j > row_i) {
758 U_cols.push_back(col_j);
759 }
760 else {
761 add_to_heap(col_j, L_cols_heap, L_cols_heaplen);
762 }
763 }
764 }
765 }
766
767 rm_heap_root(L_cols_heap, L_cols_heaplen);
768 }//end of while(L_cols_heaplen) loop
769
770
771 // Put indices and values for L into arrays and then into the L_ matrix.
772
773 // first, the original entries from the L section of A:
774 for (size_type i = 0; i < (size_type)ColIndicesA.size (); ++i) {
775 if (ColIndicesA[i] < row_i) {
776 L_tmp_idx[row_i].push_back(ColIndicesA[i]);
777 L_tmpv[row_i].push_back(cur_row[ColIndicesA[i]]);
778 pattern[ColIndicesA[i]] = UNUSED;
779 }
780 }
781
782 // next, the L entries resulting from fill:
783 for (size_type j = 0; j < L_vals_heaplen; ++j) {
784 L_tmp_idx[row_i].push_back(L_vals_heap[j]);
785 L_tmpv[row_i].push_back(cur_row[L_vals_heap[j]]);
786 pattern[L_vals_heap[j]] = UNUSED;
787 }
788
789 // L has a one on the diagonal, but we don't explicitly store
790 // it. If we don't store it, then the kernel which performs the
791 // triangular solve can assume a unit diagonal, take a short-cut
792 // and perform faster.
793
794#ifdef IFPACK2_WRITE_ILUT_FACTORS
795 for (size_type ii = 0; ii < L_tmp_idx[row_i].size (); ++ii) {
796 ofsL << row_i << " " << L_tmp_idx[row_i][ii] << " "
797 << L_tmpv[row_i][ii] << std::endl;
798 }
799#endif
800
801
802 // Pick out the diagonal element, store its reciprocal.
803 if (cur_row[row_i] == zero) {
804 std::cerr << "Ifpack2::ILUT::Compute: zero pivot encountered! "
805 << "Replacing with rownorm and continuing..."
806 << "(You may need to set the parameter "
807 << "'fact: absolute threshold'.)" << std::endl;
808 cur_row[row_i] = rownorm;
809 }
810 InvDiagU[row_i] = one / cur_row[row_i];
811
812 // Non-inverted diagonal is stored for U:
813 U_tmp_idx[row_i].push_back(row_i);
814 U_tmpv[row_i].push_back(cur_row[row_i]);
815 unorm[row_i] = scalar_mag(cur_row[row_i]);
816 pattern[row_i] = UNUSED;
817
818 // Now put indices and values for U into arrays and then into the U_ matrix.
819 // The first entry in U_cols is the diagonal, which we just handled, so we'll
820 // start our loop at j=1.
821
822 size_type U_vals_heaplen = 0;
823 for(size_type j=1; j<U_cols.size(); ++j) {
824 local_ordinal_type col = U_cols[j];
825 if (pattern[col] != ORIG) {
826 if (U_vals_heaplen < fillU) {
827 add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
828 }
829 else if (U_vals_heaplen!=0 && scalar_mag(cur_row[col]) >
830 scalar_mag(cur_row[U_vals_heap.front()])) {
831 rm_heap_root(U_vals_heap, U_vals_heaplen, vals_comp);
832 add_to_heap(col, U_vals_heap, U_vals_heaplen, vals_comp);
833 }
834 }
835 else {
836 U_tmp_idx[row_i].push_back(col);
837 U_tmpv[row_i].push_back(cur_row[col]);
838 unorm[row_i] += scalar_mag(cur_row[col]);
839 }
840 pattern[col] = UNUSED;
841 }
842
843 for(size_type j=0; j<U_vals_heaplen; ++j) {
844 U_tmp_idx[row_i].push_back(U_vals_heap[j]);
845 U_tmpv[row_i].push_back(cur_row[U_vals_heap[j]]);
846 unorm[row_i] += scalar_mag(cur_row[U_vals_heap[j]]);
847 }
848
849 unorm[row_i] /= (orig_U_len + U_vals_heaplen);
850
851#ifdef IFPACK2_WRITE_ILUT_FACTORS
852 for(int ii=0; ii<U_tmp_idx[row_i].size(); ++ii) {
853 ofsU <<row_i<< " " <<U_tmp_idx[row_i][ii]<< " "
854 <<U_tmpv[row_i][ii]<< std::endl;
855 }
856#endif
857
858 L_cols_heap.clear();
859 U_cols.clear();
860 L_vals_heap.clear();
861 U_vals_heap.clear();
862 } // end of for(row_i) loop
863
864 // Now allocate and fill the matrices
865 Array<size_t> nnzPerRow(myNumRows);
866
867 // Make sure to release the old memory for L & U prior to recomputing to
868 // avoid bloating the high-water mark.
869 L_ = Teuchos::null;
870 U_ = Teuchos::null;
871 L_solver_->setMatrix(Teuchos::null);
872 U_solver_->setMatrix(Teuchos::null);
873
874 for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
875 nnzPerRow[row_i] = L_tmp_idx[row_i].size();
876 }
877
878 L_ = rcp (new crs_matrix_type (A_local_->getRowMap(), A_local_->getColMap(),
879 nnzPerRow()));
880
881 for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
882 L_->insertLocalValues (row_i, L_tmp_idx[row_i](), L_tmpv[row_i]());
883 }
884
885 L_->fillComplete();
886
887 for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
888 nnzPerRow[row_i] = U_tmp_idx[row_i].size();
889 }
890
891 U_ = rcp (new crs_matrix_type (A_local_->getRowMap(), A_local_->getColMap(),
892 nnzPerRow()));
893
894 for (local_ordinal_type row_i = 0 ; row_i < myNumRows ; ++row_i) {
895 U_->insertLocalValues (row_i, U_tmp_idx[row_i](), U_tmpv[row_i]());
896 }
897
898 U_->fillComplete();
899
900 L_solver_->setMatrix(L_);
901 L_solver_->initialize ();
902 L_solver_->compute ();
903
904 U_solver_->setMatrix(U_);
905 U_solver_->initialize ();
906 U_solver_->compute ();
907
908 } //if (!this->useKokkosKernelsParILUT_)
909 else {
910 // Set L, U rowmaps back to original state. Par_ilut can change them, which invalidates them
911 // if compute is called again.
912 if (this->isComputed()) {
913 Kokkos::resize(L_rowmap_, L_rowmap_orig_.size());
914 Kokkos::resize(U_rowmap_, U_rowmap_orig_.size());
915 Kokkos::deep_copy(L_rowmap_, L_rowmap_orig_);
916 Kokkos::deep_copy(U_rowmap_, U_rowmap_orig_);
917 }
918
919 RCP<const crs_matrix_type> A_local_crs = Teuchos::rcp_dynamic_cast<const crs_matrix_type>(A_local_);
920 {//Make sure values in A is picked up even in case of pattern reuse
921 if(A_local_crs.is_null()) {
922 local_ordinal_type numRows = A_local_->getLocalNumRows();
923 Array<size_t> entriesPerRow(numRows);
924 for(local_ordinal_type i = 0; i < numRows; i++) {
925 entriesPerRow[i] = A_local_->getNumEntriesInLocalRow(i);
926 }
927 RCP<crs_matrix_type> A_local_crs_nc =
928 rcp (new crs_matrix_type (A_local_->getRowMap (),
929 A_local_->getColMap (),
930 entriesPerRow()));
931 // copy entries into A_local_crs
932 nonconst_local_inds_host_view_type indices("indices",A_local_->getLocalMaxNumRowEntries());
933 nonconst_values_host_view_type values("values",A_local_->getLocalMaxNumRowEntries());
934 for(local_ordinal_type i = 0; i < numRows; i++) {
935 size_t numEntries = 0;
936 A_local_->getLocalRowCopy(i, indices, values, numEntries);
937 A_local_crs_nc->insertLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()),indices.data());
938 }
939 A_local_crs_nc->fillComplete (A_local_->getDomainMap (), A_local_->getRangeMap ());
940 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
941 }
942 auto lclMtx = A_local_crs->getLocalMatrixDevice();
943 A_local_rowmap_ = lclMtx.graph.row_map;
944 A_local_entries_ = lclMtx.graph.entries;
945 A_local_values_ = lclMtx.values;
946 }
947
948 //JHU TODO Should allocation of L & U's column (aka entry) and value arrays occur here or in init()?
949 auto par_ilut_handle = KernelHandle_->get_par_ilut_handle();
950 auto nnzL = par_ilut_handle->get_nnzL();
951 static_graph_entries_t L_entries_ = static_graph_entries_t("L_entries", nnzL);
952 local_matrix_values_t L_values_ = local_matrix_values_t("L_values", nnzL);
953
954 auto nnzU = par_ilut_handle->get_nnzU();
955 static_graph_entries_t U_entries_ = static_graph_entries_t("U_entries", nnzU);
956 local_matrix_values_t U_values_ = local_matrix_values_t("U_values", nnzU);
957
958 KokkosSparse::Experimental::par_ilut_numeric(KernelHandle_.getRawPtr(),
959 A_local_rowmap_, A_local_entries_, A_local_values_,
960 L_rowmap_, L_entries_, L_values_, U_rowmap_, U_entries_, U_values_);
961
962 auto L_kokkosCrsGraph = local_graph_device_type(L_entries_, L_rowmap_);
963 auto U_kokkosCrsGraph = local_graph_device_type(U_entries_, U_rowmap_);
964
965 local_matrix_device_type L_localCrsMatrix_device;
966 L_localCrsMatrix_device = local_matrix_device_type("L_Factor_localmatrix",
967 A_local_->getLocalNumRows(),
968 L_values_,
969 L_kokkosCrsGraph);
970
971 L_ = rcp (new crs_matrix_type (L_localCrsMatrix_device,
972 A_local_crs->getRowMap(),
973 A_local_crs->getColMap(),
974 A_local_crs->getDomainMap(),
975 A_local_crs->getRangeMap(),
976 A_local_crs->getGraph()->getImporter(),
977 A_local_crs->getGraph()->getExporter()));
978
979 local_matrix_device_type U_localCrsMatrix_device;
980 U_localCrsMatrix_device = local_matrix_device_type("U_Factor_localmatrix",
981 A_local_->getLocalNumRows(),
982 U_values_,
983 U_kokkosCrsGraph);
984
985 U_ = rcp (new crs_matrix_type (U_localCrsMatrix_device,
986 A_local_crs->getRowMap(),
987 A_local_crs->getColMap(),
988 A_local_crs->getDomainMap(),
989 A_local_crs->getRangeMap(),
990 A_local_crs->getGraph()->getImporter(),
991 A_local_crs->getGraph()->getExporter()));
992
993 L_solver_->setMatrix (L_);
994 L_solver_->compute ();//NOTE: Only do compute if the pointer changed. Otherwise, do nothing
995 U_solver_->setMatrix (U_);
996 U_solver_->compute ();//NOTE: Only do compute if the pointer changed. Otherwise, do nothing
997 } //if (!this->useKokkosKernelsParILUT_) ... else ...
998
999 } // Timer scope for timing compute()
1000 ComputeTime_ += (timer.wallTime() - startTime);
1001 IsComputed_ = true;
1002 ++NumCompute_;
1003} //compute()
1004
1005
1006template <class MatrixType>
1008apply (const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
1009 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
1010 Teuchos::ETransp mode,
1011 scalar_type alpha,
1012 scalar_type beta) const
1013{
1014 using Teuchos::RCP;
1015 using Teuchos::rcp;
1016 using Teuchos::rcpFromRef;
1017
1018 TEUCHOS_TEST_FOR_EXCEPTION(
1019 ! isComputed (), std::runtime_error,
1020 "Ifpack2::ILUT::apply: You must call compute() to compute the incomplete "
1021 "factorization, before calling apply().");
1022
1023 TEUCHOS_TEST_FOR_EXCEPTION(
1024 X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
1025 "Ifpack2::ILUT::apply: X and Y must have the same number of columns. "
1026 "X has " << X.getNumVectors () << " columns, but Y has "
1027 << Y.getNumVectors () << " columns.");
1028
1029 const scalar_type one = STS::one ();
1030 const scalar_type zero = STS::zero ();
1031
1032 Teuchos::Time timer ("ILUT::apply");
1033 double startTime = timer.wallTime();
1034 { // Start timing
1035 Teuchos::TimeMonitor timeMon (timer, true);
1036
1037 if (alpha == one && beta == zero) {
1038 if (mode == Teuchos::NO_TRANS) { // Solve L (U Y) = X for Y.
1039 // Start by solving L Y = X for Y.
1040 L_solver_->apply (X, Y, mode);
1041
1042 // Solve U Y = Y.
1043 U_solver_->apply (Y, Y, mode);
1044 }
1045 else { // Solve U^P (L^P Y)) = X for Y (where P is * or T).
1046
1047 // Start by solving U^P Y = X for Y.
1048 U_solver_->apply (X, Y, mode);
1049
1050 // Solve L^P Y = Y.
1051 L_solver_->apply (Y, Y, mode);
1052 }
1053 }
1054 else { // alpha != 1 or beta != 0
1055 if (alpha == zero) {
1056 // The special case for beta == 0 ensures that if Y contains Inf
1057 // or NaN values, we replace them with 0 (following BLAS
1058 // convention), rather than multiplying them by 0 to get NaN.
1059 if (beta == zero) {
1060 Y.putScalar (zero);
1061 } else {
1062 Y.scale (beta);
1063 }
1064 } else { // alpha != zero
1065 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
1066 apply (X, Y_tmp, mode);
1067 Y.update (alpha, Y_tmp, beta);
1068 }
1069 }
1070 }//end timing
1071
1072 ++NumApply_;
1073 ApplyTime_ += (timer.wallTime() - startTime);
1074} //apply()
1075
1076
1077template <class MatrixType>
1079{
1080 std::ostringstream os;
1081
1082 // Output is a valid YAML dictionary in flow style. If you don't
1083 // like everything on a single line, you should call describe()
1084 // instead.
1085 os << "\"Ifpack2::ILUT\": {";
1086 os << "Initialized: " << (isInitialized () ? "true" : "false") << ", "
1087 << "Computed: " << (isComputed () ? "true" : "false") << ", ";
1088
1089 os << "Level-of-fill: " << getLevelOfFill() << ", "
1090 << "absolute threshold: " << getAbsoluteThreshold() << ", "
1091 << "relative threshold: " << getRelativeThreshold() << ", "
1092 << "relaxation value: " << getRelaxValue() << ", ";
1093
1094 if (A_.is_null ()) {
1095 os << "Matrix: null";
1096 }
1097 else {
1098 os << "Global matrix dimensions: ["
1099 << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"
1100 << ", Global nnz: " << A_->getGlobalNumEntries();
1101 }
1102
1103 os << "}";
1104 return os.str ();
1105}
1106
1107
1108template <class MatrixType>
1109void
1111describe (Teuchos::FancyOStream& out,
1112 const Teuchos::EVerbosityLevel verbLevel) const
1113{
1114 using Teuchos::Comm;
1115 using Teuchos::OSTab;
1116 using Teuchos::RCP;
1117 using Teuchos::TypeNameTraits;
1118 using std::endl;
1119 using Teuchos::VERB_DEFAULT;
1120 using Teuchos::VERB_NONE;
1121 using Teuchos::VERB_LOW;
1122 using Teuchos::VERB_MEDIUM;
1123 using Teuchos::VERB_HIGH;
1124 using Teuchos::VERB_EXTREME;
1125
1126 const Teuchos::EVerbosityLevel vl =
1127 (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1128 OSTab tab0 (out);
1129
1130 if (vl > VERB_NONE) {
1131 out << "\"Ifpack2::ILUT\":" << endl;
1132 OSTab tab1 (out);
1133 out << "MatrixType: " << TypeNameTraits<MatrixType>::name () << endl;
1134 if (this->getObjectLabel () != "") {
1135 out << "Label: \"" << this->getObjectLabel () << "\"" << endl;
1136 }
1137 out << "Initialized: " << (isInitialized () ? "true" : "false")
1138 << endl
1139 << "Computed: " << (isComputed () ? "true" : "false")
1140 << endl
1141 << "Level of fill: " << getLevelOfFill () << endl
1142 << "Absolute threshold: " << getAbsoluteThreshold () << endl
1143 << "Relative threshold: " << getRelativeThreshold () << endl
1144 << "Relax value: " << getRelaxValue () << endl;
1145
1146 if (isComputed () && vl >= VERB_HIGH) {
1147 const double fillFraction =
1148 (double) getGlobalNumEntries () / (double) A_->getGlobalNumEntries ();
1149 const double nnzToRows =
1150 (double) getGlobalNumEntries () / (double) U_->getGlobalNumRows ();
1151
1152 out << "Dimensions of L: [" << L_->getGlobalNumRows () << ", "
1153 << L_->getGlobalNumRows () << "]" << endl
1154 << "Dimensions of U: [" << U_->getGlobalNumRows () << ", "
1155 << U_->getGlobalNumRows () << "]" << endl
1156 << "Number of nonzeros in factors: " << getGlobalNumEntries () << endl
1157 << "Fill fraction of factors over A: " << fillFraction << endl
1158 << "Ratio of nonzeros to rows: " << nnzToRows << endl;
1159 }
1160
1161 out << "Number of initialize calls: " << getNumInitialize () << endl
1162 << "Number of compute calls: " << getNumCompute () << endl
1163 << "Number of apply calls: " << getNumApply () << endl
1164 << "Total time in seconds for initialize: " << getInitializeTime () << endl
1165 << "Total time in seconds for compute: " << getComputeTime () << endl
1166 << "Total time in seconds for apply: " << getApplyTime () << endl;
1167
1168 out << "Local matrix:" << endl;
1169 A_local_->describe (out, vl);
1170 }
1171}
1172
1173}//namespace Ifpack2
1174
1175
1176// FIXME (mfh 16 Sep 2014) We should really only use RowMatrix here!
1177// There's no need to instantiate for CrsMatrix too. All Ifpack2
1178// preconditioners can and should do dynamic casts if they need a type
1179// more specific than RowMatrix.
1180
1181#define IFPACK2_ILUT_INSTANT(S,LO,GO,N) \
1182 template class Ifpack2::ILUT< Tpetra::RowMatrix<S, LO, GO, N> >;
1183
1184#endif /* IFPACK2_ILUT_DEF_HPP */
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_ILUT_decl.hpp:87
double getLevelOfFill() const
The level of fill.
Definition Ifpack2_ILUT_decl.hpp:327
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the input matrix's communicator.
Definition Ifpack2_ILUT_def.hpp:276
double getInitializeTime() const
Returns the time spent in Initialize().
Definition Ifpack2_ILUT_def.hpp:341
void compute()
Compute factors L and U using the specified diagonal perturbation thresholds and relaxation parameter...
Definition Ifpack2_ILUT_def.hpp:553
Teuchos::RCP< const row_matrix_type > getMatrix() const
Returns a reference to the matrix to be preconditioned.
Definition Ifpack2_ILUT_def.hpp:287
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition Ifpack2_ILUT_decl.hpp:210
magnitude_type getDropTolerance() const
Gets the dropping tolerance.
Definition Ifpack2_ILUT_decl.hpp:347
global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition Ifpack2_ILUT_def.hpp:370
int getNumInitialize() const
Returns the number of calls to Initialize().
Definition Ifpack2_ILUT_def.hpp:323
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_ILUT_decl.hpp:78
Teuchos::RCP< const map_type > getDomainMap() const
Tpetra::Map representing the domain of this operator.
Definition Ifpack2_ILUT_def.hpp:294
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the ILUT preconditioner to X, resulting in Y.
Definition Ifpack2_ILUT_def.hpp:1008
int getNumApply() const
Returns the number of calls to apply().
Definition Ifpack2_ILUT_def.hpp:335
bool isComputed() const
If compute() is completed, this query returns true, otherwise it returns false.
Definition Ifpack2_ILUT_decl.hpp:225
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_ILUT_def.hpp:359
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_ILUT_decl.hpp:75
std::string description() const
Return a simple one-line description of this object.
Definition Ifpack2_ILUT_def.hpp:1078
bool hasTransposeApply() const
Whether this object's apply() method can apply the transpose (or conjugate transpose,...
Definition Ifpack2_ILUT_def.hpp:317
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
Definition Ifpack2_ILUT_def.hpp:1111
ILUT(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition Ifpack2_ILUT_def.hpp:100
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_ILUT_def.hpp:382
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Type of the Tpetra::CrsMatrix specialization that this class uses for the L and U factors.
Definition Ifpack2_ILUT_decl.hpp:109
Teuchos::RCP< const map_type > getRangeMap() const
Tpetra::Map representing the range of this operator.
Definition Ifpack2_ILUT_def.hpp:306
size_t getLocalNumEntries() const
Returns the number of nonzero entries in the local graph.
Definition Ifpack2_ILUT_def.hpp:376
void initialize()
Clear any previously computed factors, and potentially compute sparsity patterns of factors.
Definition Ifpack2_ILUT_def.hpp:455
magnitude_type getRelaxValue() const
Get the relax value.
Definition Ifpack2_ILUT_decl.hpp:342
int getNumCompute() const
Returns the number of calls to Compute().
Definition Ifpack2_ILUT_def.hpp:329
double getApplyTime() const
Returns the time spent in apply().
Definition Ifpack2_ILUT_def.hpp:353
void setParameters(const Teuchos::ParameterList &params)
Set preconditioner parameters.
Definition Ifpack2_ILUT_def.hpp:132
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition Ifpack2_ILUT_decl.hpp:337
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition Ifpack2_ILUT_decl.hpp:332
double getComputeTime() const
Returns the time spent in Compute().
Definition Ifpack2_ILUT_def.hpp:347
Access only local rows and columns of a sparse matrix.
Definition Ifpack2_LocalFilter_decl.hpp:131
"Preconditioner" that solves local sparse triangular systems.
Definition Ifpack2_LocalSparseTriangularSolver_decl.hpp:55
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41
void add_to_heap(const Ordinal &idx, Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition Ifpack2_Heap.hpp:37
void rm_heap_root(Teuchos::Array< Ordinal > &heap, SizeType &heap_len)
Definition Ifpack2_Heap.hpp:59