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>
27 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
32 , computedFlag_(false)
40 , params_(Params::getDefaults()) {}
42 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
46 return mat_->getDomainMap();
49 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
53 return mat_->getRangeMap();
56 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
58 apply(
const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
59 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Y,
63 const std::string timerName(
"Ifpack2::FastILU::apply");
70 if (!isInitialized() || !isComputed()) {
71 throw std::runtime_error(std::string(
"Called ") + getName() +
"::apply() without first calling initialize() and/or compute().");
73 if (X.getNumVectors() != Y.getNumVectors()) {
74 throw std::invalid_argument(getName() +
"::apply: X and Y have different numbers of vectors (pass X and Y with exactly matching dimensions)");
76 if (X.getLocalLength() != Y.getLocalLength()) {
77 throw std::invalid_argument(getName() +
"::apply: X and Y have different lengths (pass X and Y with exactly matching dimensions)");
81 int nvecs = X.getNumVectors();
82 auto nrowsX = X.getLocalLength();
83 auto nrowsY = Y.getLocalLength();
85 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
86 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
90 applyLocalPrec(x1d, y1d);
93 auto x2d = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
94 auto y2d = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
95 for (
int i = 0; i < nvecs; i++) {
96 auto xColView1d = Kokkos::subview(x2d, Kokkos::ALL(), i);
97 auto yColView1d = Kokkos::subview(y2d, Kokkos::ALL(), i);
98 ImplScalarArray x1d(const_cast<ImplScalar*>(xColView1d.data()), nrowsX);
99 ImplScalarArray y1d(const_cast<ImplScalar*>(yColView1d.data()), nrowsY);
101 applyLocalPrec(x1d, y1d);
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 params_ = Params(List, getName());
113 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 return params_.blockCrs && params_.blockCrsSize > 1;
119 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
122 const std::string timerName(
"Ifpack2::FastILU::initialize");
129 if (mat_.is_null()) {
130 throw std::runtime_error(std::string(
"Called ") + getName() +
"::initialize() but matrix was null (call setMatrix() with a non-null matrix first)");
134 auto crs_matrix = Ifpack2::Details::getCrsMatrix(this->mat_);
136 if (params_.fillBlocks) {
138 auto crs_matrix_block_filled = Tpetra::fillLogicalBlocks(*crs_matrix, params_.blockCrsSize);
139 auto bcrs_matrix = Tpetra::convertToBlockCrsMatrix(*crs_matrix_block_filled, params_.blockCrsSize);
143 auto bcrs_matrix = Tpetra::convertToBlockCrsMatrix(*crs_matrix, params_.blockCrsSize);
148 Kokkos::Timer copyTimer;
149 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getStructure(mat_.get(), localRowPtrsHost_, localRowPtrs_, localColInds_);
150 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
151 crsCopyTime_ = copyTimer.seconds();
153 if (params_.use_metis) {
154 assert(!params_.blockCrs);
155 const std::string timerNameMetis(
"Ifpack2::FastILU::Metis");
161 #ifdef HAVE_IFPACK2_METIS
162 idx_t nrows = localRowPtrsHost_.size() - 1;
165 metis_perm_ = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing(
"metis_perm"), nrows);
166 metis_iperm_ = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing(
"metis_iperm"), nrows);
169 auto localColIndsHost_ = Kokkos::create_mirror_view(localColInds_);
170 Kokkos::deep_copy(localColIndsHost_, localColInds_);
173 idx_t nnz = localColIndsHost_.size();
174 MetisArrayHost metis_rowptr;
175 MetisArrayHost metis_colidx;
177 bool metis_symmetrize =
true;
178 if (metis_symmetrize) {
180 using OrdinalArrayMirror =
typename OrdinalArray::host_mirror_type;
181 KokkosKernels::Impl::symmetrize_graph_symbolic_hashmap<
182 OrdinalArrayHost, OrdinalArrayMirror, MetisArrayHost, MetisArrayHost, Kokkos::HostSpace::execution_space>(nrows, localRowPtrsHost_, localColIndsHost_, metis_rowptr, metis_colidx);
185 idx_t old_nnz = nnz = 0;
186 for (idx_t i = 0; i < nrows; i++) {
187 for (LocalOrdinal k = old_nnz; k < metis_rowptr(i + 1); k++) {
188 if (metis_colidx(k) != i) {
189 metis_colidx(nnz) = metis_colidx(k);
193 old_nnz = metis_rowptr(i + 1);
194 metis_rowptr(i + 1) = nnz;
198 metis_rowptr = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing(
"metis_rowptr"), nrows + 1);
199 metis_colidx = MetisArrayHost(Kokkos::ViewAllocateWithoutInitializing(
"metis_colidx"), nnz);
202 for (idx_t i = 0; i < nrows; i++) {
203 for (LocalOrdinal k = localRowPtrsHost_(i); k < localRowPtrsHost_(i + 1); k++) {
204 if (localColIndsHost_(k) != i) {
205 metis_colidx(nnz) = localColIndsHost_(k);
209 metis_rowptr(i + 1) = nnz;
214 int info = METIS_NodeND(&nrows, metis_rowptr.data(), metis_colidx.data(),
215 NULL, NULL, metis_perm_.data(), metis_iperm_.data());
216 if (METIS_OK != info) {
217 throw std::runtime_error(std::string(
"METIS_NodeND returned info = " + info));
221 throw std::runtime_error(std::string(
"TPL METIS is not enabled"));
230 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
240 throw std::runtime_error(getName() +
": initialize() must be called before compute()");
243 const std::string timerName(
"Ifpack2::FastILU::compute");
251 Kokkos::Timer copyTimer;
252 CrsArrayReader<Scalar, ImplScalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
253 crsCopyTime_ += copyTimer.seconds();
255 computedFlag_ =
true;
259 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
262 return computedFlag_;
265 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
272 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
278 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
284 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
290 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
296 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
308 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
314 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
318 throw std::runtime_error(std::string(
"Preconditioner type Ifpack2::Details::") + getName() +
" doesn't support checkLocalILU().");
321 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 throw std::runtime_error(std::string(
"Preconditioner type Ifpack2::Details::") + getName() +
" doesn't support checkLocalIC().");
328 template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
330 std::ostringstream os;
332 os <<
"\"Ifpack2::Details::" << getName() <<
"\": {";
333 os <<
"Initialized: " << (isInitialized() ?
"true" :
"false") <<
", ";
334 os <<
"Computed: " << (isComputed() ?
"true" :
"false") <<
", ";
335 os <<
"Sweeps: " << getSweeps() <<
", ";
336 os <<
"Triangular solve type: " << getSpTrsvType() <<
", ";
337 if (getSpTrsvType() ==
"Fast") {
338 os <<
"# of triangular solve iterations: " << getNTrisol() <<
", ";
340 if (mat_.is_null()) {
341 os <<
"Matrix: null";
343 os <<
"Global matrix dimensions: [" << mat_->getGlobalNumRows() <<
", " << mat_->getGlobalNumCols() <<
"]";
344 os <<
", Global nnz: " << mat_->getGlobalNumEntries();
349 template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
353 throw std::invalid_argument(std::string(
"Ifpack2::Details::") + getName() +
"::setMatrix() called with a null matrix. Pass a non-null matrix.");
356 if (mat_.get() != A.
get()) {
359 computedFlag_ =
false;
363 template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
369 p.sptrsv_algo = FastILU::SpTRSV::Fast;
380 p.fillBlocks =
false;
384 template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
385 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
387 *
this = getDefaults();
392 #define TYPE_ERROR(name, correctTypeName) \
393 { throw std::invalid_argument(precType + "::setParameters(): parameter \"" + name + "\" has the wrong type (must be " + correctTypeName + ")"); }
394 #define CHECK_VALUE(param, member, cond, msg) \
397 throw std::invalid_argument(precType + "::setParameters(): parameter \"" + param + "\" has value " + std::to_string(member) + " but " + msg); \
403 if (pL.
isType<
bool>(
"metis"))
404 use_metis = pL.
get<
bool>(
"metis");
406 TYPE_ERROR(
"metis",
"bool");
410 if (pL.
isType<
int>(
"sweeps")) {
411 nFact = pL.
get<
int>(
"sweeps");
412 CHECK_VALUE(
"sweeps", nFact, nFact < 1,
"must have a value of at least 1");
414 TYPE_ERROR(
"sweeps",
"int");
416 std::string sptrsv_type =
"Fast";
418 sptrsv_type = pL.
get<std::string>(
"triangular solve type");
420 if (sptrsv_type ==
"Standard Host") {
421 sptrsv_algo = FastILU::SpTRSV::StandardHost;
422 }
else if (sptrsv_type ==
"Standard") {
423 sptrsv_algo = FastILU::SpTRSV::Standard;
427 if (pL.
isParameter(
"triangular solve iterations")) {
428 if (pL.
isType<
int>(
"triangular solve iterations")) {
429 nTrisol = pL.
get<
int>(
"triangular solve iterations");
430 CHECK_VALUE(
"triangular solve iterations", nTrisol, nTrisol < 1,
"must have a value of at least 1");
432 TYPE_ERROR(
"triangular solve iterations",
"int");
436 if (pL.
isType<
int>(
"level")) {
437 level = pL.
get<
int>(
"level");
438 }
else if (pL.
isType<
double>(
"level")) {
441 double dval = pL.
get<
double>(
"level");
443 double fpart = modf(dval, &ipart);
445 CHECK_VALUE(
"level", level, fpart != 0,
"must be an integral value");
447 TYPE_ERROR(
"level",
"int");
449 CHECK_VALUE(
"level", level, level < 0,
"must be nonnegative");
452 if (pL.
isType<
double>(
"damping factor"))
453 omega = pL.
get<
double>(
"damping factor");
455 TYPE_ERROR(
"damping factor",
"double");
458 if (pL.
isType<
double>(
"shift"))
459 shift = pL.
get<
double>(
"shift");
461 TYPE_ERROR(
"shift",
"double");
465 if (pL.
isType<
bool>(
"guess"))
466 guessFlag = pL.
get<
bool>(
"guess");
468 TYPE_ERROR(
"guess",
"bool");
472 if (pL.
isType<
int>(
"block size for ILU")) {
473 blockSizeILU = pL.
get<
int>(
"block size for ILU");
474 CHECK_VALUE(
"block size for ILU", blockSizeILU, blockSizeILU < 1,
"must have a value of at least 1");
476 TYPE_ERROR(
"block size for ILU",
"int");
480 if (pL.
isType<
int>(
"block size for SpTRSV"))
481 blockSize = pL.
get<
int>(
"block size for SpTRSV");
483 TYPE_ERROR(
"block size for SpTRSV",
"int");
487 if (pL.
isType<
bool>(
"block crs"))
488 blockCrs = pL.
get<
bool>(
"block crs");
490 TYPE_ERROR(
"block crs",
"bool");
494 if (pL.
isType<
int>(
"block crs block size"))
495 blockCrsSize = pL.
get<
int>(
"block crs block size");
497 TYPE_ERROR(
"block crs block size",
"int");
501 if (pL.
isType<
bool>(
"fill blocks for input"))
502 blockCrsSize = pL.
get<
bool>(
"fill blocks for input");
504 TYPE_ERROR(
"fill blocks for input",
"bool");
511 #define IFPACK2_DETAILS_FASTILU_BASE_INSTANT(S, L, G, N) \
512 template class Ifpack2::Details::FastILU_Base<S, L, G, N>;
int getNumCompute() const
Get the number of times compute() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:280
virtual void checkLocalIC() const
Verify and print debug information about the underlying IC preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:323
double getComputeTime() const
Get the time spent in the last compute() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:298
T & get(const std::string &name, T def_value)
double getCopyTime() const
Get the time spent deep copying local 3-array CRS out of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:310
double getInitializeTime() const
Get the time spent in the last initialize() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:292
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Get the range map of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:52
double getApplyTime() const
Get the time spent in the last apply() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:304
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Get the domain map of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:45
bool isParameter(const std::string &name) const
Teuchos::RCP< const TRowMatrix > getMatrix() const
Get the current matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:268
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:316
int getNumApply() const
Get the number of times apply() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:286
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:36
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:108
std::string description() const
Return a brief description of the preconditioner, in YAML format.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:329
void compute()
Compute the preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:238
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:29
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:261
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:232
void initialize()
Initialize the preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:121
bool isType(const std::string &name) const
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition: Ifpack2_Details_FastILU_Base_def.hpp:351
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:59
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:58
int getNumInitialize() const
Get the number of times initialize() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:274
Kokkos::View< LocalOrdinal *, execution_space >::HostMirror OrdinalArrayHost
Array of LocalOrdinal on host.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:57
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...