Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_Filu_def.hpp
1/*@HEADER
2// ***********************************************************************
3//
4// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5// Copyright (2009) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ***********************************************************************
40//@HEADER
41*/
42
44
45#ifndef __IFPACK2_FILU_DEF_HPP__
46#define __IFPACK2_FILU_DEF_HPP__
47
48#include "Ifpack2_Details_Filu_decl.hpp"
50#include "Ifpack2_Details_getCrsMatrix.hpp"
51#include <Kokkos_Timer.hpp>
52#include <shylu_fastilu.hpp>
53
54namespace Ifpack2
55{
56namespace Details
57{
58
59template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
61Filu(Teuchos::RCP<const TRowMatrix> A) :
62 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
63
64template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
66getSweeps() const
67{
68 return localPrec_->getNFact();
69}
70
71template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
73getSpTrsvType() const
74{
75 return localPrec_->getSpTrsvType();
76}
77
78template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
80getNTrisol() const
81{
82 return localPrec_->getNTrisol();
83}
84
85template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
87checkLocalILU() const
88{
89 localPrec_->checkILU();
90}
91
92template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
94checkLocalIC() const
95{
96 localPrec_->checkIC();
97}
98
99template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
102{
103 auto nRows = this->mat_->getLocalNumRows();
104 auto& p = this->params_;
105 auto matCrs = Ifpack2::Details::getCrsMatrix(this->mat_);
106
107 bool skipSortMatrix = !matCrs.is_null() && matCrs->getCrsGraph()->isSorted() &&
108 !p.use_metis;
109 localPrec_ = Teuchos::rcp(new LocalFILU(skipSortMatrix, this->localRowPtrs_, this->localColInds_, this->localValues_, nRows, p.sptrsv_algo,
110 p.nFact, p.nTrisol, p.level, p.omega, p.shift, p.guessFlag ? 1 : 0, p.blockSizeILU, p.blockSize));
111 #ifdef HAVE_IFPACK2_METIS
112 if (p.use_metis) {
113 localPrec_->setMetisPerm(this->metis_perm_, this->metis_iperm_);
114 }
115 #endif
116 localPrec_->initialize();
117 this->initTime_ = localPrec_->getInitializeTime();
118}
119
120template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
123{
124 //update values in local prec (until compute(), values aren't needed)
125 localPrec_->setValues(this->localValues_);
126 localPrec_->compute();
127 this->computeTime_ = localPrec_->getComputeTime();
128}
129
130template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
132applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
133{
134 localPrec_->apply(x, y);
135 //since this may be applied to multiple vectors, add to applyTime_ instead of setting it
136 this->applyTime_ += localPrec_->getApplyTime();
137}
138
139template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
141getName() const
142{
143 return "Filu";
144}
145
146#define IFPACK2_DETAILS_FILU_INSTANT(S, L, G, N) \
147template class Ifpack2::Details::Filu<S, L, G, N>;
148
149} //namespace Details
150} //namespace Ifpack2
151
152#endif
153
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:74
void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only)
Definition: Ifpack2_Details_Filu_def.hpp:132
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:66
void checkLocalILU() const
Verify and print debug info about the internal ILU preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:87
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition: Ifpack2_Details_Filu_def.hpp:141
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:80
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:122
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:94
Filu(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Filu_def.hpp:61
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition: Ifpack2_Details_Filu_def.hpp:101
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition: Ifpack2_Details_Filu_def.hpp:73
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:74