Ifpack2 Templated Preconditioning Package Version 1.0
Ifpack2_Details_FastILU_Base_def.hpp
Go to the documentation of this file.
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_FASTILU_BASE_DEF_HPP__
46#define __IFPACK2_FASTILU_BASE_DEF_HPP__
47
49#include <impl/Kokkos_Timer.hpp>
50#include <stdexcept>
51#include "Teuchos_TimeMonitor.hpp"
52
53
54namespace Ifpack2
55{
56namespace Details
57{
58
59template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61FastILU_Base(Teuchos::RCP<const TRowMatrix> A) :
62 mat_(A),
63 initFlag_(false),
64 computedFlag_(false),
65 nInit_(0),
66 nComputed_(0),
67 nApply_(0),
68 initTime_(0.0),
69 computeTime_(0.0),
70 applyTime_(0.0),
71 crsCopyTime_(0.0),
72 params_(Params::getDefaults()) {}
73
74template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
77getDomainMap () const
78{
79 return mat_->getDomainMap();
80}
81
82template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
85getRangeMap () const
86{
87 return mat_->getRangeMap();
88}
89
90template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92apply (const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
93 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
94 Teuchos::ETransp mode,
95 Scalar alpha,
96 Scalar beta) const
97{
98 const std::string timerName ("Ifpack2::FastILU::apply");
99 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
100 if (timer.is_null ()) {
101 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
102 }
103 Teuchos::TimeMonitor timeMon (*timer);
104
105 if(!isInitialized() || !isComputed())
106 {
107 throw std::runtime_error(std::string("Called ") + getName() + "::apply() without first calling initialize() and/or compute().");
108 }
109 if(X.getNumVectors() != Y.getNumVectors())
110 {
111 throw std::invalid_argument(getName() + "::apply: X and Y have different numbers of vectors (pass X and Y with exactly matching dimensions)");
112 }
113 if(X.getLocalLength() != Y.getLocalLength())
114 {
115 throw std::invalid_argument(getName() + "::apply: X and Y have different lengths (pass X and Y with exactly matching dimensions)");
116 }
117 //zero out applyTime_ now, because the calls to applyLocalPrec() will add to it
118 applyTime_ = 0;
119 int nvecs = X.getNumVectors();
120 if(nvecs == 1)
121 {
122 auto x2d = X.template getLocalView<execution_space>(Tpetra::Access::ReadWrite);
123 auto y2d = Y.template getLocalView<execution_space>(Tpetra::Access::ReadWrite);
124 auto x1d = Kokkos::subview(x2d, Kokkos::ALL(), 0);
125 auto y1d = Kokkos::subview(y2d, Kokkos::ALL(), 0);
126 applyLocalPrec(x1d, y1d);
127 }
128 else
129 {
130 //Solve each vector one at a time (until FastILU supports multiple RHS)
131 for(int i = 0; i < nvecs; i++)
132 {
133 auto Xcol = X.getVector(i);
134 auto Ycol = Y.getVector(i);
135 auto xColView2d = Xcol->template getLocalView<execution_space>(Tpetra::Access::ReadWrite);
136 auto yColView2d = Ycol->template getLocalView<execution_space>(Tpetra::Access::ReadWrite);
137 ScalarArray xColView1d = Kokkos::subview(xColView2d, Kokkos::ALL(), 0);
138 ScalarArray yColView1d = Kokkos::subview(yColView2d, Kokkos::ALL(), 0);
139 applyLocalPrec(xColView1d, yColView1d);
140 }
141 }
142}
143
144template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146setParameters (const Teuchos::ParameterList& List)
147{
148 //Params constructor does all parameter validation, and sets default values
149 params_ = Params(List, getName());
150}
151
152template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
155{
156
157 const std::string timerName ("Ifpack2::FastILU::initialize");
158 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
159 if (timer.is_null ()) {
160 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
161 }
162
163 if(mat_.is_null())
164 {
165 throw std::runtime_error(std::string("Called ") + getName() + "::initialize() but matrix was null (call setMatrix() with a non-null matrix first)");
166 }
167 Kokkos::Impl::Timer copyTimer;
168 CrsArrayReader<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getStructure(mat_.get(), localRowPtrsHost_, localRowPtrs_, localColInds_);
169 crsCopyTime_ = copyTimer.seconds();
170 initLocalPrec(); //note: initLocalPrec updates initTime
171 initFlag_ = true;
172 nInit_++;
173}
174
175template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177isInitialized() const
178{
179 return initFlag_;
180}
181
182template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
184compute()
185{
186 if(!initFlag_)
187 {
188 throw std::runtime_error(getName() + ": initialize() must be called before compute()");
189 }
190
191 const std::string timerName ("Ifpack2::FastILU::compute");
192 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter (timerName);
193 if (timer.is_null ()) {
194 timer = Teuchos::TimeMonitor::getNewCounter (timerName);
195 }
196
197
198 //get copy of values array from matrix
199 Kokkos::Impl::Timer copyTimer;
200 CrsArrayReader<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getValues(mat_.get(), localValues_, localRowPtrsHost_);
201 crsCopyTime_ += copyTimer.seconds(); //add to the time spent getting rowptrs/colinds
202 computeLocalPrec(); //this updates computeTime_
203 computedFlag_ = true;
204 nComputed_++;
205}
206
207template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
209isComputed() const
210{
211 return computedFlag_;
212}
213
214template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
215Teuchos::RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
217getMatrix() const
218{
219 return mat_;
220}
221
222template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224getNumInitialize() const
225{
226 return nInit_;
227}
228
229template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231getNumCompute() const
232{
233 return nComputed_;
234}
235
236template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
238getNumApply() const
239{
240 return nApply_;
241}
242
243template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
245getInitializeTime() const
246{
247 return initTime_;
248}
249
250template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
252getComputeTime() const
253{
254 return computeTime_;
255}
256
257template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
259getApplyTime() const
260{
261 return applyTime_;
262}
263
264template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266getCopyTime() const
267{
268 return crsCopyTime_;
269}
270
271template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
273checkLocalILU() const
274{
275 //if the underlying type of this doesn't implement checkLocalILU, it's an illegal operation
276 throw std::runtime_error(std::string("Preconditioner type Ifpack2::Details::") + getName() + " doesn't support checkLocalILU().");
277}
278
279template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281checkLocalIC() const
282{
283 //if the underlying type of this doesn't implement checkLocalIC, it's an illegal operation
284 throw std::runtime_error(std::string("Preconditioner type Ifpack2::Details::") + getName() + " doesn't support checkLocalIC().");
285}
286
287template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
289{
290 std::ostringstream os;
291 //Output is a YAML dictionary
292 os << "\"Ifpack2::Details::" << getName() << "\": {";
293 os << "Initialized: " << (isInitialized() ? "true" : "false") << ", ";
294 os << "Computed: " << (isComputed() ? "true" : "false") << ", ";
295 os << "Sweeps: " << getSweeps() << ", ";
296 os << "# of triangular solve iterations: " << getNTrisol() << ", ";
297 if(mat_.is_null())
298 {
299 os << "Matrix: null";
300 }
301 else
302 {
303 os << "Global matrix dimensions: [" << mat_->getGlobalNumRows() << ", " << mat_->getGlobalNumCols() << "]";
304 os << ", Global nnz: " << mat_->getGlobalNumEntries();
305 }
306 return os.str();
307}
308
309template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
311setMatrix(const Teuchos::RCP<const TRowMatrix>& A)
312{
313 if(A.is_null())
314 {
315 throw std::invalid_argument(std::string("Ifpack2::Details::") + getName() + "::setMatrix() called with a null matrix. Pass a non-null matrix.");
316 }
317 typedef Tpetra::RowGraph<LocalOrdinal, GlobalOrdinal, Node> RGraph;
318 Teuchos::RCP<const RGraph> aGraph; //graph of A
319 Teuchos::RCP<const RGraph> matGraph; //graph of current mat_
320 try
321 {
322 aGraph = A->getGraph();
323 }
324 catch(...)
325 {
326 aGraph = Teuchos::null;
327 }
328 if(!mat_.is_null())
329 {
330 try
331 {
332 matGraph = mat_->getGraph();
333 }
334 catch(...)
335 {
336 matGraph = Teuchos::null;
337 }
338 }
339 //bmk note: this modeled after RILUK::setMatrix
340 if(mat_.get() != A.get())
341 {
342 mat_ = A;
343 if(matGraph.is_null() || (matGraph.getRawPtr() != aGraph.getRawPtr()))
344 {
345 //must assume that matrix's graph changed, so need to copy the structure again in initialize()
346 initFlag_ = false;
347 }
348 computedFlag_ = false;
349 }
350}
351
352template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
356{
357 Params p;
358 p.nFact = 5;
359 p.nTrisol = 1;
360 p.level = 0;
361 p.omega = 0.5;
362 p.shift = 0;
363 p.guessFlag = true;
364 p.blockSize = 1;
365 return p;
366}
367
368template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
369FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
370Params::Params(const Teuchos::ParameterList& pL, std::string precType)
371{
372 *this = getDefaults();
373 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude;
374 //For each parameter, check that if the parameter exists, it has the right type
375 //Then get the value and sanity check it
376 //If the parameter doesn't exist, leave it as default (from Params::getDefaults())
377 //"sweeps" aka nFact
378 #define TYPE_ERROR(name, correctTypeName) {throw std::invalid_argument(precType + "::setParameters(): parameter \"" + name + "\" has the wrong type (must be " + correctTypeName + ")");}
379 #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);}}
380 if(pL.isParameter("sweeps"))
381 {
382 if(pL.isType<int>("sweeps"))
383 {
384 nFact = pL.get<int>("sweeps");
385 CHECK_VALUE("sweeps", nFact, nFact < 1, "must have a value of at least 1");
386 }
387 else
388 TYPE_ERROR("sweeps", "int");
389 }
390 //"triangular solve iterations" aka nTrisol
391 if(pL.isParameter("triangular solve iterations"))
392 {
393 if(pL.isType<int>("triangular solve iterations"))
394 {
395 nTrisol = pL.get<int>("triangular solve iterations");
396 CHECK_VALUE("triangular solve iterations", nTrisol, nTrisol < 1, "must have a value of at least 1");
397 }
398 else
399 TYPE_ERROR("triangular solve iterations", "int");
400 }
401 //"level"
402 if(pL.isParameter("level"))
403 {
404 if(pL.isType<int>("level"))
405 {
406 level = pL.get<int>("level");
407 }
408 else if(pL.isType<double>("level"))
409 {
410 //Level can be read as double (like in ILUT), but must be an exact integer
411 //Any integer used for level-of-fill can be represented exactly in double (so use exact compare)
412 double dval = pL.get<double>("level");
413 double ipart;
414 double fpart = modf(dval, &ipart);
415 level = ipart;
416 CHECK_VALUE("level", level, fpart != 0, "must be an integral value");
417 }
418 else
419 {
420 TYPE_ERROR("level", "int");
421 }
422 CHECK_VALUE("level", level, level < 0, "must be nonnegative");
423 }
424 //"damping factor" aka omega -- try both double and magnitude as type
425 if(pL.isParameter("damping factor"))
426 {
427 if(pL.isType<double>("damping factor"))
428 omega = pL.get<double>("damping factor");
429 else if(pL.isType<magnitude>("damping factor"))
430 omega = pL.get<magnitude>("damping factor");
431 else
432 TYPE_ERROR("damping factor", "double or magnitude_type");
433 CHECK_VALUE("damping factor", omega, omega <= 0 || omega > 1, "must be in the range (0, 1]");
434 }
435 //"shift" -- also try both double and magnitude
436 if(pL.isParameter("shift"))
437 {
438 if(pL.isType<double>("shift"))
439 shift = pL.get<double>("shift");
440 else if(pL.isType<magnitude>("shift"))
441 shift = pL.get<magnitude>("shift");
442 else
443 TYPE_ERROR("shift", "double or magnitude_type");
444 //no hard requirements for shift value so don't
445 }
446 //"guess" aka guessFlag
447 if(pL.isParameter("guess"))
448 {
449 if(pL.isType<bool>("guess"))
450 guessFlag = pL.get<bool>("guess");
451 else
452 TYPE_ERROR("guess", "bool");
453 }
454 //"block size" aka blkSz
455 if(pL.isParameter("block size"))
456 {
457 if(pL.isType<int>("block size"))
458 blockSize = pL.get<int>("block size");
459 else
460 TYPE_ERROR("block size", "int");
461 }
462 #undef CHECK_VALUE
463 #undef TYPE_ERROR
464}
465
466#define IFPACK2_DETAILS_FASTILU_BASE_INSTANT(S, L, G, N) \
467template class Ifpack2::Details::FastILU_Base<S, L, G, N>;
468
469} //namespace Details
470} //namespace Ifpack2
471
472#endif
473
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:69
double getInitializeTime() const
Get the time spent in the last initialize() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:245
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:273
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition: Ifpack2_Details_FastILU_Base_def.hpp:311
void compute()
Compute the preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:184
double getComputeTime() const
Get the time spent in the last compute() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:252
double getCopyTime() const
Get the time spent deep copying local 3-array CRS out of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:266
int getNumApply() const
Get the number of times apply() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:238
void initialize()
Initialize the preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:154
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:177
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:61
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:146
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:92
double getApplyTime() const
Get the time spent in the last apply() call.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:259
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Get the domain map of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:77
virtual void checkLocalIC() const
Verify and print debug information about the underlying IC preconditioner.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:281
Kokkos::View< Scalar *, execution_space > ScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:88
int getNumInitialize() const
Get the number of times initialize() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:224
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Get the range map of the matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:85
Teuchos::RCP< const TRowMatrix > getMatrix() const
Get the current matrix.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:217
std::string description() const
Return a brief description of the preconditioner, in YAML format.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:288
int getNumCompute() const
Get the number of times compute() was called.
Definition: Ifpack2_Details_FastILU_Base_def.hpp:231
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:209
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:73