12#ifndef __IFPACK2_FASTILU_BASE_DEF_HPP__
13#define __IFPACK2_FASTILU_BASE_DEF_HPP__
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>
29template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42 params_(Params::getDefaults()) {}
44template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
45Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
49 return mat_->getDomainMap();
52template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
53Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
57 return mat_->getRangeMap();
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,
68 const std::string
timerName (
"Ifpack2::FastILU::apply");
69 Teuchos::RCP<Teuchos::Time>
timer = Teuchos::TimeMonitor::lookupCounter (
timerName);
70 if (
timer.is_null ()) {
75 if(!isInitialized() || !isComputed())
77 throw std::runtime_error(std::string(
"Called ") + getName() +
"::apply() without first calling initialize() and/or compute().");
79 if(X.getNumVectors() != Y.getNumVectors())
81 throw std::invalid_argument(getName() +
"::apply: X and Y have different numbers of vectors (pass X and Y with exactly matching dimensions)");
83 if(X.getLocalLength() != Y.getLocalLength())
85 throw std::invalid_argument(getName() +
"::apply: X and Y have different lengths (pass X and Y with exactly matching dimensions)");
89 int nvecs = X.getNumVectors();
90 auto nrowsX = X.getLocalLength();
91 auto nrowsY = Y.getLocalLength();
94 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
95 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
104 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
105 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
118template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
123 params_ = Params(
List, getName());
126template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 return params_.blockCrs && params_.blockCrsSize > 1;
133template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 const std::string
timerName (
"Ifpack2::FastILU::initialize");
138 Teuchos::RCP<Teuchos::Time>
timer = Teuchos::TimeMonitor::lookupCounter (
timerName);
139 if (
timer.is_null ()) {
146 throw std::runtime_error(std::string(
"Called ") + getName() +
"::initialize() but matrix was null (call setMatrix() with a non-null matrix first)");
150 auto crs_matrix = Ifpack2::Details::getCrsMatrix(this->mat_);
152 if (params_.fillBlocks) {
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_);
170 if (params_.use_metis)
172 assert(!params_.blockCrs);
179 #ifdef HAVE_IFPACK2_METIS
199 KokkosKernels::Impl::symmetrize_graph_symbolic_hashmap<
222 for (LocalOrdinal
k = localRowPtrsHost_(
i);
k < localRowPtrsHost_(
i+1);
k++) {
236 throw std::runtime_error(std::string(
"METIS_NodeND returned info = " +
info));
240 throw std::runtime_error(std::string(
"TPL METIS is not enabled"));
249template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
256template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
262 throw std::runtime_error(getName() +
": initialize() must be called before compute()");
265 const std::string
timerName (
"Ifpack2::FastILU::compute");
266 Teuchos::RCP<Teuchos::Time>
timer = Teuchos::TimeMonitor::lookupCounter (
timerName);
267 if (
timer.is_null ()) {
274 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
277 computedFlag_ =
true;
281template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
285 return computedFlag_;
288template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
289Teuchos::RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
296template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
303template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
310template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
317template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
324template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
331template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
338template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
345template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
350 throw std::runtime_error(std::string(
"Preconditioner type Ifpack2::Details::") + getName() +
" doesn't support checkLocalILU().");
353template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
358 throw std::runtime_error(std::string(
"Preconditioner type Ifpack2::Details::") + getName() +
" doesn't support checkLocalIC().");
361template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
364 std::ostringstream
os;
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() <<
", ";
376 os <<
"Matrix: null";
380 os <<
"Global matrix dimensions: [" << mat_->getGlobalNumRows() <<
", " << mat_->getGlobalNumCols() <<
"]";
381 os <<
", Global nnz: " << mat_->getGlobalNumEntries();
386template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
388setMatrix(
const Teuchos::RCP<const TRowMatrix>& A)
392 throw std::invalid_argument(std::string(
"Ifpack2::Details::") + getName() +
"::setMatrix() called with a null matrix. Pass a non-null matrix.");
395 if(mat_.get() != A.get())
399 computedFlag_ =
false;
403template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
410 p.sptrsv_algo = FastILU::SpTRSV::Fast;
421 p.fillBlocks =
false;
425template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
426FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
427Params::Params(
const Teuchos::ParameterList& pL, std::string precType)
429 *
this = getDefaults();
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);}}
438 if(pL.isParameter(
"metis"))
440 if(pL.isType<
bool>(
"metis"))
441 use_metis = pL.get<
bool>(
"metis");
443 TYPE_ERROR(
"metis",
"bool");
446 if(pL.isParameter(
"sweeps"))
448 if(pL.isType<
int>(
"sweeps"))
450 nFact = pL.get<
int>(
"sweeps");
451 CHECK_VALUE(
"sweeps", nFact, nFact < 1,
"must have a value of at least 1");
454 TYPE_ERROR(
"sweeps",
"int");
456 std::string sptrsv_type =
"Fast";
457 if(pL.isParameter(
"triangular solve type")) {
458 sptrsv_type = pL.get<std::string> (
"triangular solve type");
460 if (sptrsv_type ==
"Standard Host") {
461 sptrsv_algo = FastILU::SpTRSV::StandardHost;
462 }
else if (sptrsv_type ==
"Standard") {
463 sptrsv_algo = FastILU::SpTRSV::Standard;
467 if(pL.isParameter(
"triangular solve iterations"))
469 if(pL.isType<
int>(
"triangular solve iterations"))
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");
475 TYPE_ERROR(
"triangular solve iterations",
"int");
478 if(pL.isParameter(
"level"))
480 if(pL.isType<
int>(
"level"))
482 level = pL.get<
int>(
"level");
484 else if(pL.isType<
double>(
"level"))
488 double dval = pL.get<
double>(
"level");
490 double fpart = modf(dval, &ipart);
492 CHECK_VALUE(
"level", level, fpart != 0,
"must be an integral value");
496 TYPE_ERROR(
"level",
"int");
498 CHECK_VALUE(
"level", level, level < 0,
"must be nonnegative");
500 if(pL.isParameter(
"damping factor"))
502 if(pL.isType<
double>(
"damping factor"))
503 omega = pL.get<
double>(
"damping factor");
505 TYPE_ERROR(
"damping factor",
"double");
507 if(pL.isParameter(
"shift"))
509 if(pL.isType<
double>(
"shift"))
510 shift = pL.get<
double>(
"shift");
512 TYPE_ERROR(
"shift",
"double");
515 if(pL.isParameter(
"guess"))
517 if(pL.isType<
bool>(
"guess"))
518 guessFlag = pL.get<
bool>(
"guess");
520 TYPE_ERROR(
"guess",
"bool");
523 if(pL.isParameter(
"block size for ILU"))
525 if(pL.isType<
int>(
"block size for ILU"))
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");
531 TYPE_ERROR(
"block size for ILU",
"int");
534 if(pL.isParameter(
"block size for SpTRSV"))
536 if(pL.isType<
int>(
"block size for SpTRSV"))
537 blockSize = pL.get<
int>(
"block size for SpTRSV");
539 TYPE_ERROR(
"block size for SpTRSV",
"int");
542 if(pL.isParameter(
"block crs"))
544 if(pL.isType<
bool>(
"block crs"))
545 blockCrs = pL.get<
bool>(
"block crs");
547 TYPE_ERROR(
"block crs",
"bool");
550 if(pL.isParameter(
"block crs block size"))
552 if(pL.isType<
int>(
"block crs block size"))
553 blockCrsSize = pL.get<
int>(
"block crs block size");
555 TYPE_ERROR(
"block crs block size",
"int");
558 if(pL.isParameter(
"fill blocks for input"))
560 if(pL.isType<
bool>(
"fill blocks for input"))
561 blockCrsSize = pL.get<
bool>(
"fill blocks for input");
563 TYPE_ERROR(
"fill blocks for input",
"bool");
570#define IFPACK2_DETAILS_FASTILU_BASE_INSTANT(S, L, G, N) \
571template class Ifpack2::Details::FastILU_Base<S, L, G, N>;
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
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition Ifpack2_Details_FastILU_Base_def.hpp:388
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
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_FastILU_Base_def.hpp:31
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
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
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
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:77
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41