Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_Fildl_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_FILDL_DEF_HPP__
46#define __IFPACK2_FILDL_DEF_HPP__
47
48#include "Ifpack2_Details_Fildl_decl.hpp"
50#include <Kokkos_Timer.hpp>
51#include <shylu_fastildl.hpp>
52
53namespace Ifpack2
54{
55namespace Details
56{
57
58template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
60Fildl(Teuchos::RCP<const TRowMatrix> A) :
61 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
62
63template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
65getSweeps() const
66{
67 return localPrec_->getNFact();
68}
69
70template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
72getSpTrsvType() const
73{
74 return localPrec_->getSpTrsvType();
75}
76
77template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
79getNTrisol() const
80{
81 return localPrec_->getNTrisol();
82}
83
84template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
86checkLocalIC() const
87{
88 localPrec_->checkIC();
89}
90
91template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
94{
95 auto nRows = this->mat_->getLocalNumRows();
96 auto& p = this->params_;
97 localPrec_ = Teuchos::rcp(new LocalFILDL(this->localRowPtrs_, this->localColInds_, this->localValues_, nRows, (p.sptrsv_algo != FastILU::SpTRSV::Fast),
98 p.nFact, p.nTrisol, p.level, p.omega, p.shift, p.guessFlag ? 1 : 0, p.blockSizeILU, p.blockSize));
99 localPrec_->initialize();
100 this->initTime_ = localPrec_->getInitializeTime();
101}
102
103template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
106{
107 //update values in local prec (until compute(), values aren't needed)
108 localPrec_->setValues(this->localValues_);
109 localPrec_->compute();
110 this->computeTime_ = localPrec_->getComputeTime();
111}
112
113template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
115applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
116{
117 localPrec_->apply(x, y);
118 //since this may be applied to multiple vectors, add to applyTime_ instead of setting it
119 this->applyTime_ += localPrec_->getApplyTime();
120}
121
122template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
124getName() const
125{
126 return "Fildl";
127}
128
129#define IFPACK2_DETAILS_FILDL_INSTANT(S, L, G, N) \
130template class Ifpack2::Details::Fildl<S, L, G, N>;
131
132} //namespace Details
133} //namespace Ifpack2
134
135#endif
136
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_Fildl_def.hpp:115
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:105
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:65
Fildl(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_Fildl_def.hpp:60
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition Ifpack2_Details_Fildl_def.hpp:124
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Fildl_def.hpp:72
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:79
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:86
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Fildl_def.hpp:93
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74