Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 
54 namespace Ifpack2
55 {
56 namespace Details
57 {
58 
59 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
62  FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
63 
64 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
66 getSweeps() const
67 {
68  return localPrec_->getNFact();
69 }
70 
71 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
74 {
75  return localPrec_->getSpTrsvType();
76 }
77 
78 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
80 getNTrisol() const
81 {
82  return localPrec_->getNTrisol();
83 }
84 
85 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
88 {
89  localPrec_->checkILU();
90 }
91 
92 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
94 checkLocalIC() const
95 {
96  localPrec_->checkIC();
97 }
98 
99 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
102 {
103  auto nRows = this->mat_->getLocalNumRows();
104  auto& p = this->params_;
105  auto matCrs = Ifpack2::Details::getCrsMatrix(this->mat_);
106 
107  if (p.blockCrsSize > 1 && !BlockCrsEnabled) {
108  throw std::runtime_error("Must use prec type FAST_ILU_B if you want blockCrs support");
109  }
110 
111  bool skipSortMatrix = !matCrs.is_null() && matCrs->getCrsGraph()->isSorted() &&
112  !p.use_metis;
113  localPrec_ =
114  Teuchos::rcp(new LocalFILU(skipSortMatrix, this->localRowPtrs_, this->localColInds_, this->localValues_, nRows, p.sptrsv_algo,
115  p.nFact, p.nTrisol, p.level, p.omega, p.shift, p.guessFlag ? 1 : 0, p.blockSizeILU, p.blockSize,
116  p.blockCrsSize));
117 
118  #ifdef HAVE_IFPACK2_METIS
119  if (p.use_metis) {
120  localPrec_->setMetisPerm(this->metis_perm_, this->metis_iperm_);
121  }
122  #endif
123 
124  localPrec_->initialize();
125  this->initTime_ = localPrec_->getInitializeTime();
126 }
127 
128 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
131 {
132  //update values in local prec (until compute(), values aren't needed)
133  localPrec_->setValues(this->localValues_);
134  localPrec_->compute();
135  this->computeTime_ = localPrec_->getComputeTime();
136 }
137 
138 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
141 {
142  localPrec_->apply(x, y);
143 
144  //since this may be applied to multiple vectors, add to applyTime_ instead of setting it
145  this->applyTime_ += localPrec_->getApplyTime();
146 }
147 
148 template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
150 getName() const
151 {
152  return "Filu";
153 }
154 
155 #define IFPACK2_DETAILS_FILU_INSTANT(S, L, G, N) \
156 template class Ifpack2::Details::Filu<S, L, G, N,false>; \
157 template class Ifpack2::Details::Filu<S, L, G, N,true>;
158 
159 } //namespace Details
160 } //namespace Ifpack2
161 
162 #endif
163 
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition: Ifpack2_Details_Filu_def.hpp:73
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:66
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:130
std::string getName() const
Get the name of the underlying preconditioner (&quot;Filu&quot;, &quot;Fildl&quot; or &quot;Fic&quot;)
Definition: Ifpack2_Details_Filu_def.hpp:150
void checkLocalILU() const
Verify and print debug info about the internal ILU preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:87
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_-&gt;ini...
Definition: Ifpack2_Details_Filu_def.hpp:101
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:80
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:71
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:95
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:140
Filu(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Filu_def.hpp:61
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:94
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...