Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_FastILU_Base_def.hpp
Go to the documentation of this file.
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
11
12#ifndef __IFPACK2_FASTILU_BASE_DEF_HPP__
13#define __IFPACK2_FASTILU_BASE_DEF_HPP__
14
16#include "Tpetra_BlockCrsMatrix.hpp"
17#include "Tpetra_BlockCrsMatrix_Helpers.hpp"
18#include "Ifpack2_Details_getCrsMatrix.hpp"
19#include <KokkosKernels_Utils.hpp>
20#include <Kokkos_Timer.hpp>
21#include <Teuchos_TimeMonitor.hpp>
22#include <stdexcept>
23
24namespace Ifpack2
25{
26namespace Details
27{
28
29template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31FastILU_Base(Teuchos::RCP<const TRowMatrix> A) :
32 mat_(A),
33 initFlag_(false),
34 computedFlag_(false),
35 nInit_(0),
36 nComputed_(0),
37 nApply_(0),
38 initTime_(0.0),
39 computeTime_(0.0),
40 applyTime_(0.0),
41 crsCopyTime_(0.0),
42 params_(Params::getDefaults()) {}
43
44template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
47getDomainMap () const
48{
49 return mat_->getDomainMap();
50}
51
52template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
53Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
55getRangeMap () const
56{
57 return mat_->getRangeMap();
58}
59
60template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62apply (const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
63 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
64 Teuchos::ETransp mode,
65 Scalar alpha,
66 Scalar beta) const
67{
68 const std::string timerName ("Ifpack2::FastILU::apply");
69 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
70 if (timer.is_null ()) {
71 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
72 }
73 Teuchos::TimeMonitor timeMon (*timer);
74
75 if(!isInitialized() || !isComputed())
76 {
77 throw std::runtime_error(std::string("Called ") + getName() + "::apply() without first calling initialize() and/or compute().");
78 }
79 if(X.getNumVectors() != Y.getNumVectors())
80 {
81 throw std::invalid_argument(getName() + "::apply: X and Y have different numbers of vectors (pass X and Y with exactly matching dimensions)");
82 }
83 if(X.getLocalLength() != Y.getLocalLength())
84 {
85 throw std::invalid_argument(getName() + "::apply: X and Y have different lengths (pass X and Y with exactly matching dimensions)");
86 }
87 //zero out applyTime_ now, because the calls to applyLocalPrec() will add to it
88 applyTime_ = 0;
89 int nvecs = X.getNumVectors();
90 auto nrowsX = X.getLocalLength();
91 auto nrowsY = Y.getLocalLength();
92 if(nvecs == 1)
93 {
94 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
95 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
96 ImplScalarArray x1d (const_cast<ImplScalar*>(x2d.data()), nrowsX);
97 ImplScalarArray y1d (const_cast<ImplScalar*>(y2d.data()), nrowsY);
98
99 applyLocalPrec(x1d, y1d);
100 }
101 else
102 {
103 //Solve each vector one at a time (until FastILU supports multiple RHS)
104 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
105 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
106 for(int i = 0; i < nvecs; i++)
107 {
108 auto xColView1d = Kokkos::subview(x2d, Kokkos::ALL(), i);
109 auto yColView1d = Kokkos::subview(y2d, Kokkos::ALL(), i);
110 ImplScalarArray x1d (const_cast<ImplScalar*>(xColView1d.data()), nrowsX);
111 ImplScalarArray y1d (const_cast<ImplScalar*>(yColView1d.data()), nrowsY);
112
113 applyLocalPrec(x1d, y1d);
114 }
115 }
116}
117
118template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120setParameters (const Teuchos::ParameterList& List)
121{
122 //Params constructor does all parameter validation, and sets default values
123 params_ = Params(List, getName());
124}
125
126template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
128isBlockCrs() const
129{
130 return params_.blockCrs && params_.blockCrsSize > 1;
131}
132
133template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
136{
137 const std::string timerName ("Ifpack2::FastILU::initialize");
138 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
139 if (timer.is_null ()) {
140 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
141 }
142 Teuchos::TimeMonitor timeMon (*timer);
143
144 if(mat_.is_null())
145 {
146 throw std::runtime_error(std::string("Called ") + getName() + "::initialize() but matrix was null (call setMatrix() with a non-null matrix first)");
147 }
148
149 if (isBlockCrs()) {
150 auto crs_matrix = Ifpack2::Details::getCrsMatrix(this->mat_);
151
152 if (params_.fillBlocks) {
153 // Create new TCrsMatrix with the new filled data and conver to Bcrs
154 auto crs_matrix_block_filled = Tpetra::fillLogicalBlocks(*crs_matrix, params_.blockCrsSize);
155 auto bcrs_matrix = Tpetra::convertToBlockCrsMatrix(*crs_matrix_block_filled, params_.blockCrsSize);
156 mat_ = bcrs_matrix;
157 }
158 else {
159 // Assume input is already filled, just convert to Bcrs
160 auto bcrs_matrix = Tpetra::convertToBlockCrsMatrix(*crs_matrix, params_.blockCrsSize);
161 mat_ = bcrs_matrix;
162 }
163 }
164
165 Kokkos::Timer copyTimer;
166 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getStructure(mat_.get(), localRowPtrsHost_, localRowPtrs_, localColInds_);
167 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
168 crsCopyTime_ = copyTimer.seconds();
169
170 if (params_.use_metis)
171 {
172 assert(!params_.blockCrs);
173 const std::string timerNameMetis ("Ifpack2::FastILU::Metis");
174 Teuchos::RCP<Teuchos::Time> timerMetis = Teuchos::TimeMonitor::lookupCounter (timerNameMetis);
175 if (timerMetis.is_null ()) {
176 timerMetis = Teuchos::TimeMonitor::getNewCounter (timerNameMetis);
177 }
178 Teuchos::TimeMonitor timeMonMetis (*timerMetis);
179 #ifdef HAVE_IFPACK2_METIS
180 idx_t nrows = localRowPtrsHost_.size() - 1;
181 if (nrows > 0) {
182 // reorder will convert both graph and perm/iperm to the internal METIS integer type
183 metis_perm_ = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing("metis_perm"), nrows);
184 metis_iperm_ = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing("metis_iperm"), nrows);
185
186 // copy ColInds to host
187 auto localColIndsHost_ = Kokkos::create_mirror_view(localColInds_);
188 Kokkos::deep_copy(localColIndsHost_, localColInds_);
189
190 // prepare for calling metis
191 idx_t nnz = localColIndsHost_.size();
192 MetisArrayHost metis_rowptr;
193 MetisArrayHost metis_colidx;
194
195 bool metis_symmetrize = true;
196 if (metis_symmetrize) {
197 // symmetrize
198 using OrdinalArrayMirror = typename OrdinalArray::host_mirror_type;
199 KokkosKernels::Impl::symmetrize_graph_symbolic_hashmap<
200 OrdinalArrayHost, OrdinalArrayMirror, MetisArrayHost, MetisArrayHost, Kokkos::HostSpace::execution_space>
201 (nrows, localRowPtrsHost_, localColIndsHost_, metis_rowptr, metis_colidx);
202
203 // remove diagonals
204 idx_t old_nnz = nnz = 0;
205 for (idx_t i = 0; i < nrows; i++) {
206 for (LocalOrdinal k = old_nnz; k < metis_rowptr(i+1); k++) {
207 if (metis_colidx(k) != i) {
208 metis_colidx(nnz) = metis_colidx(k);
209 nnz++;
210 }
211 }
212 old_nnz = metis_rowptr(i+1);
213 metis_rowptr(i+1) = nnz;
214 }
215 } else {
216 // copy and remove diagonals
217 metis_rowptr = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing("metis_rowptr"), nrows+1);
218 metis_colidx = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing("metis_colidx"), nnz);
219 nnz = 0;
220 metis_rowptr(0) = 0;
221 for (idx_t i = 0; i < nrows; i++) {
222 for (LocalOrdinal k = localRowPtrsHost_(i); k < localRowPtrsHost_(i+1); k++) {
223 if (localColIndsHost_(k) != i) {
224 metis_colidx(nnz) = localColIndsHost_(k);
225 nnz++;
226 }
227 }
228 metis_rowptr(i+1) = nnz;
229 }
230 }
231
232 // call metis
233 int info = METIS_NodeND(&nrows, metis_rowptr.data(), metis_colidx.data(),
234 NULL, NULL, metis_perm_.data(), metis_iperm_.data());
235 if (METIS_OK != info) {
236 throw std::runtime_error(std::string("METIS_NodeND returned info = " + info));
237 }
238 }
239 #else
240 throw std::runtime_error(std::string("TPL METIS is not enabled"));
241 #endif
242 }
243
244 initLocalPrec(); //note: initLocalPrec updates initTime
245 initFlag_ = true;
246 nInit_++;
247}
248
249template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
251isInitialized() const
252{
253 return initFlag_;
254}
255
256template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258compute()
259{
260 if(!initFlag_)
261 {
262 throw std::runtime_error(getName() + ": initialize() must be called before compute()");
263 }
264
265 const std::string timerName ("Ifpack2::FastILU::compute");
266 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
267 if (timer.is_null ()) {
268 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
269 }
270 Teuchos::TimeMonitor timeMon (*timer);
271
272 //get copy of values array from matrix
273 Kokkos::Timer copyTimer;
274 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
275 crsCopyTime_ += copyTimer.seconds(); //add to the time spent getting rowptrs/colinds
276 computeLocalPrec(); //this updates computeTime_
277 computedFlag_ = true;
278 nComputed_++;
279}
280
281template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
283isComputed() const
284{
285 return computedFlag_;
286}
287
288template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289Teuchos::RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
295
296template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
302
303template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305getNumCompute() const
306{
307 return nComputed_;
308}
309
310template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
316
317template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
319getInitializeTime() const
320{
321 return initTime_;
322}
323
324template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326getComputeTime() const
327{
328 return computeTime_;
329}
330
331template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
333getApplyTime() const
334{
335 return applyTime_;
336}
337
338template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
340getCopyTime() const
341{
342 return crsCopyTime_;
343}
344
345template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
347checkLocalILU() const
348{
349 //if the underlying type of this doesn't implement checkLocalILU, it's an illegal operation
350 throw std::runtime_error(std::string("Preconditioner type Ifpack2::Details::") + getName() + " doesn't support checkLocalILU().");
351}
352
353template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
355checkLocalIC() const
356{
357 //if the underlying type of this doesn't implement checkLocalIC, it's an illegal operation
358 throw std::runtime_error(std::string("Preconditioner type Ifpack2::Details::") + getName() + " doesn't support checkLocalIC().");
359}
360
361template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
363{
364 std::ostringstream os;
365 //Output is a YAML dictionary
366 os << "\"Ifpack2::Details::" << getName() << "\": {";
367 os << "Initialized: " << (isInitialized() ? "true" : "false") << ", ";
368 os << "Computed: " << (isComputed() ? "true" : "false") << ", ";
369 os << "Sweeps: " << getSweeps() << ", ";
370 os << "Triangular solve type: " << getSpTrsvType() << ", ";
371 if (getSpTrsvType() == "Fast") {
372 os << "# of triangular solve iterations: " << getNTrisol() << ", ";
373 }
374 if(mat_.is_null())
375 {
376 os << "Matrix: null";
377 }
378 else
379 {
380 os << "Global matrix dimensions: [" << mat_->getGlobalNumRows() << ", " << mat_->getGlobalNumCols() << "]";
381 os << ", Global nnz: " << mat_->getGlobalNumEntries();
382 }
383 return os.str();
384}
385
386template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
388setMatrix(const Teuchos::RCP<const TRowMatrix>& A)
389{
390 if(A.is_null())
391 {
392 throw std::invalid_argument(std::string("Ifpack2::Details::") + getName() + "::setMatrix() called with a null matrix. Pass a non-null matrix.");
393 }
394 //bmk note: this modeled after RILUK::setMatrix
395 if(mat_.get() != A.get())
396 {
397 mat_ = A;
398 initFlag_ = false;
399 computedFlag_ = false;
400 }
401}
402
403template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
404typename FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Params
407{
408 Params p;
409 p.use_metis = false;
410 p.sptrsv_algo = FastILU::SpTRSV::Fast;
411 p.nFact = 5; // # of sweeps for computing fastILU
412 p.nTrisol = 5; // # of sweeps for applying fastSpTRSV
413 p.level = 0; // level of ILU
414 p.omega = 1.0; // damping factor for fastILU
415 p.shift = 0;
416 p.guessFlag = true;
417 p.blockSizeILU = 1; // # of nonzeros / thread, for fastILU
418 p.blockSize = 1; // # of rows / thread, for SpTRSV
419 p.blockCrs = false; // whether to use block CRS for fastILU
420 p.blockCrsSize = 1; // block size for block CRS
421 p.fillBlocks = false; // whether input matrix needs to be filled
422 return p;
423}
424
425template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
426FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
427Params::Params(const Teuchos::ParameterList& pL, std::string precType)
428{
429 *this = getDefaults();
430 //For each parameter, check that if the parameter exists, it has the right type
431 //Then get the value and sanity check it
432 //If the parameter doesn't exist, leave it as default (from Params::getDefaults())
433 //"sweeps" aka nFact
434 #define TYPE_ERROR(name, correctTypeName) {throw std::invalid_argument(precType + "::setParameters(): parameter \"" + name + "\" has the wrong type (must be " + correctTypeName + ")");}
435 #define CHECK_VALUE(param, member, cond, msg) {if(cond) {throw std::invalid_argument(precType + "::setParameters(): parameter \"" + param + "\" has value " + std::to_string(member) + " but " + msg);}}
436
437 //metis
438 if(pL.isParameter("metis"))
439 {
440 if(pL.isType<bool>("metis"))
441 use_metis = pL.get<bool>("metis");
442 else
443 TYPE_ERROR("metis", "bool");
444 }
445
446 if(pL.isParameter("sweeps"))
447 {
448 if(pL.isType<int>("sweeps"))
449 {
450 nFact = pL.get<int>("sweeps");
451 CHECK_VALUE("sweeps", nFact, nFact < 1, "must have a value of at least 1");
452 }
453 else
454 TYPE_ERROR("sweeps", "int");
455 }
456 std::string sptrsv_type = "Fast";
457 if(pL.isParameter("triangular solve type")) {
458 sptrsv_type = pL.get<std::string> ("triangular solve type");
459 }
460 if (sptrsv_type == "Standard Host") {
461 sptrsv_algo = FastILU::SpTRSV::StandardHost;
462 } else if (sptrsv_type == "Standard") {
463 sptrsv_algo = FastILU::SpTRSV::Standard;
464 }
465
466 //"triangular solve iterations" aka nTrisol
467 if(pL.isParameter("triangular solve iterations"))
468 {
469 if(pL.isType<int>("triangular solve iterations"))
470 {
471 nTrisol = pL.get<int>("triangular solve iterations");
472 CHECK_VALUE("triangular solve iterations", nTrisol, nTrisol < 1, "must have a value of at least 1");
473 }
474 else
475 TYPE_ERROR("triangular solve iterations", "int");
476 }
477 //"level"
478 if(pL.isParameter("level"))
479 {
480 if(pL.isType<int>("level"))
481 {
482 level = pL.get<int>("level");
483 }
484 else if(pL.isType<double>("level"))
485 {
486 //Level can be read as double (like in ILUT), but must be an exact integer
487 //Any integer used for level-of-fill can be represented exactly in double (so use exact compare)
488 double dval = pL.get<double>("level");
489 double ipart;
490 double fpart = modf(dval, &ipart);
491 level = ipart;
492 CHECK_VALUE("level", level, fpart != 0, "must be an integral value");
493 }
494 else
495 {
496 TYPE_ERROR("level", "int");
497 }
498 CHECK_VALUE("level", level, level < 0, "must be nonnegative");
499 }
500 if(pL.isParameter("damping factor"))
501 {
502 if(pL.isType<double>("damping factor"))
503 omega = pL.get<double>("damping factor");
504 else
505 TYPE_ERROR("damping factor", "double");
506 }
507 if(pL.isParameter("shift"))
508 {
509 if(pL.isType<double>("shift"))
510 shift = pL.get<double>("shift");
511 else
512 TYPE_ERROR("shift", "double");
513 }
514 //"guess" aka guessFlag
515 if(pL.isParameter("guess"))
516 {
517 if(pL.isType<bool>("guess"))
518 guessFlag = pL.get<bool>("guess");
519 else
520 TYPE_ERROR("guess", "bool");
521 }
522 //"block size" aka blkSz
523 if(pL.isParameter("block size for ILU"))
524 {
525 if(pL.isType<int>("block size for ILU"))
526 {
527 blockSizeILU = pL.get<int>("block size for ILU");
528 CHECK_VALUE("block size for ILU", blockSizeILU, blockSizeILU < 1, "must have a value of at least 1");
529 }
530 else
531 TYPE_ERROR("block size for ILU", "int");
532 }
533 //"block size" aka blkSz
534 if(pL.isParameter("block size for SpTRSV"))
535 {
536 if(pL.isType<int>("block size for SpTRSV"))
537 blockSize = pL.get<int>("block size for SpTRSV");
538 else
539 TYPE_ERROR("block size for SpTRSV", "int");
540 }
541 //"block crs" aka blockCrs
542 if(pL.isParameter("block crs"))
543 {
544 if(pL.isType<bool>("block crs"))
545 blockCrs = pL.get<bool>("block crs");
546 else
547 TYPE_ERROR("block crs", "bool");
548 }
549 //"block crs block size" aka blockCrsSize
550 if(pL.isParameter("block crs block size"))
551 {
552 if(pL.isType<int>("block crs block size"))
553 blockCrsSize = pL.get<int>("block crs block size");
554 else
555 TYPE_ERROR("block crs block size", "int");
556 }
557 //"fill blocks for input" aka fillBlocks
558 if(pL.isParameter("fill blocks for input"))
559 {
560 if(pL.isType<bool>("fill blocks for input"))
561 blockCrsSize = pL.get<bool>("fill blocks for input");
562 else
563 TYPE_ERROR("fill blocks for input", "bool");
564 }
565
566 #undef CHECK_VALUE
567 #undef TYPE_ERROR
568}
569
570#define IFPACK2_DETAILS_FASTILU_BASE_INSTANT(S, L, G, N) \
571template class Ifpack2::Details::FastILU_Base<S, L, G, N>;
572
573} //namespace Details
574} //namespace Ifpack2
575
576#endif
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic).
Definition Ifpack2_Details_FastILU_Base_decl.hpp:41
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type ImplScalar
Kokkos scalar type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:48
double getInitializeTime() const
Get the time spent in the last initialize() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:319
virtual void checkLocalILU() const
Verify and print debug information about the underlying ILU preconditioner (only supported if this is...
Definition Ifpack2_Details_FastILU_Base_def.hpp:347
virtual std::string getSpTrsvType() const =0
Get the name of triangular solve algorithm.
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition Ifpack2_Details_FastILU_Base_def.hpp:388
virtual void computeLocalPrec()=0
Get values array from the matrix and then call compute() on the underlying preconditioner.
void compute()
Compute the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:258
double getComputeTime() const
Get the time spent in the last compute() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:326
double getCopyTime() const
Get the time spent deep copying local 3-array CRS out of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:340
int getNumApply() const
Get the number of times apply() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:312
void initialize()
Initialize the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:135
bool isInitialized() const
Whether initialize() has been called since the last time the matrix's structure was changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:251
virtual std::string getName() const =0
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic").
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_FastILU_Base_def.hpp:31
virtual int getSweeps() const =0
Get the "sweeps" parameter.
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided.
Definition Ifpack2_Details_FastILU_Base_def.hpp:120
void apply(const TMultiVec &X, TMultiVec &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Apply the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:62
double getApplyTime() const
Get the time spent in the last apply() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:333
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Get the domain map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:47
virtual void checkLocalIC() const
Verify and print debug information about the underlying IC preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:355
virtual void initLocalPrec()=0
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
int getNumInitialize() const
Get the number of times initialize() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:298
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Get the range map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:55
Teuchos::RCP< const TRowMatrix > getMatrix() const
Get the current matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:291
Kokkos::View< LocalOrdinal *, execution_space >::HostMirror OrdinalArrayHost
Array of LocalOrdinal on host.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:60
std::string description() const
Return a brief description of the preconditioner, in YAML format.
Definition Ifpack2_Details_FastILU_Base_def.hpp:362
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:62
virtual int getNTrisol() const =0
Get the "triangular solve iterations" parameter.
int getNumCompute() const
Get the number of times compute() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:305
bool isComputed() const
Whether compute() has been called since the last time the matrix's values or structure were changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:283
virtual void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const =0
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only).
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41