Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_Chebyshev_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_DETAILS_CHEBYSHEV_DEF_HPP
11#define IFPACK2_DETAILS_CHEBYSHEV_DEF_HPP
12
19
23// #include "Ifpack2_Details_ScaledDampedResidual.hpp"
24#include "Ifpack2_Details_ChebyshevKernel.hpp"
25#include "Kokkos_ArithTraits.hpp"
26#include "Teuchos_FancyOStream.hpp"
27#include "Teuchos_oblackholestream.hpp"
28#include "Tpetra_Details_residual.hpp"
29#include "Teuchos_LAPACK.hpp"
30#include "Ifpack2_Details_LapackSupportsScalar.hpp"
31#include <cmath>
32#include <iostream>
33
34namespace Ifpack2 {
35namespace Details {
36
37namespace { // (anonymous)
38
39// We use this text a lot in error messages.
40const char computeBeforeApplyReminder[] =
41 "This means one of the following:\n"
42 " - you have not yet called compute() on this instance, or \n"
43 " - you didn't call compute() after calling setParameters().\n\n"
44 "After creating an Ifpack2::Chebyshev instance,\n"
45 "you must _always_ call compute() at least once before calling apply().\n"
46 "After calling compute() once, you do not need to call it again,\n"
47 "unless the matrix has changed or you have changed parameters\n"
48 "(by calling setParameters()).";
49
50} // namespace (anonymous)
51
52// ReciprocalThreshold stuff below needs to be in a namspace visible outside
53// of this file
54template<class XV, class SizeType = typename XV::size_type>
55struct V_ReciprocalThresholdSelfFunctor
56{
57 typedef typename XV::execution_space execution_space;
58 typedef typename XV::non_const_value_type value_type;
59 typedef SizeType size_type;
60 typedef Kokkos::ArithTraits<value_type> KAT;
61 typedef typename KAT::mag_type mag_type;
62
63 XV X_;
64 const value_type minVal_;
65 const mag_type minValMag_;
66
67 V_ReciprocalThresholdSelfFunctor (const XV& X,
68 const value_type& min_val) :
69 X_ (X),
70 minVal_ (min_val),
71 minValMag_ (KAT::abs (min_val))
72 {}
73
74 KOKKOS_INLINE_FUNCTION
75 void operator() (const size_type& i) const
76 {
77 const mag_type X_i_abs = KAT::abs (X_[i]);
78
79 if (X_i_abs < minValMag_) {
80 X_[i] = minVal_;
81 }
82 else {
83 X_[i] = KAT::one () / X_[i];
84 }
85 }
86};
87
88template<class XV, class SizeType = typename XV::size_type>
89struct LocalReciprocalThreshold {
90 static void
91 compute (const XV& X,
92 const typename XV::non_const_value_type& minVal)
93 {
94 typedef typename XV::execution_space execution_space;
95 Kokkos::RangePolicy<execution_space, SizeType> policy (0, X.extent (0));
96 V_ReciprocalThresholdSelfFunctor<XV, SizeType> op (X, minVal);
97 Kokkos::parallel_for (policy, op);
98 }
99};
100
101template <class TpetraVectorType,
102 const bool classic = TpetraVectorType::node_type::classic>
103struct GlobalReciprocalThreshold {};
104
105template <class TpetraVectorType>
106struct GlobalReciprocalThreshold<TpetraVectorType, true> {
107 static void
108 compute (TpetraVectorType& V,
109 const typename TpetraVectorType::scalar_type& min_val)
110 {
111 typedef typename TpetraVectorType::scalar_type scalar_type;
112 typedef typename TpetraVectorType::mag_type mag_type;
113 typedef Kokkos::ArithTraits<scalar_type> STS;
115 const scalar_type ONE = STS::one ();
116 const mag_type min_val_abs = STS::abs (min_val);
117
118 Teuchos::ArrayRCP<scalar_type> V_0 = V.getDataNonConst (0);
119 const size_t lclNumRows = V.getLocalLength ();
120
121 for (size_t i = 0; i < lclNumRows; ++i) {
122 const scalar_type V_0i = V_0[i];
123 if (STS::abs (V_0i) < min_val_abs) {
124 V_0[i] = min_val;
125 } else {
126 V_0[i] = ONE / V_0i;
127 }
128 }
129 }
130};
131
132template <class TpetraVectorType>
133struct GlobalReciprocalThreshold<TpetraVectorType, false> {
134 static void
135 compute (TpetraVectorType& X,
136 const typename TpetraVectorType::scalar_type& minVal)
137 {
138 typedef typename TpetraVectorType::impl_scalar_type value_type;
139
140 const value_type minValS = static_cast<value_type> (minVal);
141 auto X_0 = Kokkos::subview (X.getLocalViewDevice (Tpetra::Access::ReadWrite),
142 Kokkos::ALL (), 0);
143 LocalReciprocalThreshold<decltype (X_0) >::compute (X_0, minValS);
144 }
145};
146
147// Utility function for inverting diagonal with threshold.
148template <typename S, typename L, typename G, typename N>
149void
150reciprocal_threshold (Tpetra::Vector<S,L,G,N>& V, const S& minVal)
151{
152 GlobalReciprocalThreshold<Tpetra::Vector<S,L,G,N> >::compute (V, minVal);
153}
154
155
156template<class ScalarType, const bool lapackSupportsScalarType = LapackSupportsScalar<ScalarType>::value>
157struct LapackHelper{
158 static
159 ScalarType
160 tri_diag_spectral_radius(Teuchos::ArrayRCP<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> diag,
161 Teuchos::ArrayRCP<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> offdiag) {
162 throw std::runtime_error("LAPACK does not support the scalar type.");
163 }
164
165};
166
167template<class V>
168void
169computeInitialGuessForCG (const V& diagonal, V& x) {
170 using device_type = typename V::node_type::device_type;
171 using range_policy = Kokkos::RangePolicy<typename device_type::execution_space>;
172
173 // Initial randomization of the vector
174 x.randomize();
175
176
177 // Zero the stuff that where the diagonal is equal to one. These are assumed to
178 // correspond to OAZ rows in the matrix.
179 size_t N = x.getLocalLength();
180 auto d_view = diagonal.template getLocalView<device_type>(Tpetra::Access::ReadOnly);
181 auto x_view = x.template getLocalView<device_type>(Tpetra::Access::ReadWrite);
182
183 auto ONE = Teuchos::ScalarTraits<typename V::impl_scalar_type>::one();
184 auto ZERO = Teuchos::ScalarTraits<typename V::impl_scalar_type>::zero();
185
186 Kokkos::parallel_for("computeInitialGuessforCG::zero_bcs", range_policy(0,N), KOKKOS_LAMBDA(const size_t & i) {
187 if(d_view(i,0) == ONE)
188 x_view(i,0) = ZERO;
189 });
190}
191
192
193
194
195
196template<class ScalarType>
197struct LapackHelper<ScalarType,true> {
198 static
199 ScalarType
200 tri_diag_spectral_radius(Teuchos::ArrayRCP<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> diag,
201 Teuchos::ArrayRCP<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> offdiag) {
202 using STS = Teuchos::ScalarTraits<ScalarType>;
203 using MagnitudeType = typename STS::magnitudeType;
204 int info = 0;
205 const int N = diag.size ();
206 ScalarType scalar_dummy;
207 std::vector<MagnitudeType> mag_dummy(4*N);
208 char char_N = 'N';
209
210 // lambdaMin = one;
211 ScalarType lambdaMax = STS::one();
212 if( N > 2 ) {
213 Teuchos::LAPACK<int,ScalarType> lapack;
214 lapack.PTEQR (char_N, N, diag.getRawPtr (), offdiag.getRawPtr (),
215 &scalar_dummy, 1, &mag_dummy[0], &info);
216 TEUCHOS_TEST_FOR_EXCEPTION
217 (info < 0, std::logic_error, "Ifpack2::Details::LapackHelper::tri_diag_spectral_radius:"
218 "LAPACK's _PTEQR failed with info = "
219 << info << " < 0. This suggests there might be a bug in the way Ifpack2 "
220 "is calling LAPACK. Please report this to the Ifpack2 developers.");
221 // lambdaMin = Teuchos::as<ScalarType> (diag[N-1]);
222 lambdaMax = Teuchos::as<ScalarType> (diag[0]);
223 }
224 return lambdaMax;
225 }
226};
227
228
229template<class ScalarType, class MV>
230void Chebyshev<ScalarType, MV>::checkInputMatrix () const
231{
232 TEUCHOS_TEST_FOR_EXCEPTION(
233 ! A_.is_null () && A_->getGlobalNumRows () != A_->getGlobalNumCols (),
234 std::invalid_argument,
235 "Ifpack2::Chebyshev: The input matrix A must be square. "
236 "A has " << A_->getGlobalNumRows () << " rows and "
237 << A_->getGlobalNumCols () << " columns.");
238
239 // In debug mode, test that the domain and range Maps of the matrix
240 // are the same.
241 if (debug_ && ! A_.is_null ()) {
242 Teuchos::RCP<const map_type> domainMap = A_->getDomainMap ();
243 Teuchos::RCP<const map_type> rangeMap = A_->getRangeMap ();
244
245 // isSameAs is a collective, but if the two pointers are the same,
246 // isSameAs will assume that they are the same on all processes, and
247 // return true without an all-reduce.
248 TEUCHOS_TEST_FOR_EXCEPTION(
249 ! domainMap->isSameAs (*rangeMap), std::invalid_argument,
250 "Ifpack2::Chebyshev: The domain Map and range Map of the matrix must be "
251 "the same (in the sense of isSameAs())." << std::endl << "We only check "
252 "for this in debug mode.");
253 }
254}
255
256
257template<class ScalarType, class MV>
258void
261{
262 // mfh 12 Aug 2016: The if statement avoids an "unreachable
263 // statement" warning for the checkInputMatrix() call, when
264 // STS::isComplex is false.
265 if (STS::isComplex) {
266 TEUCHOS_TEST_FOR_EXCEPTION
267 (true, std::logic_error, "Ifpack2::Chebyshev: This class' implementation "
268 "of Chebyshev iteration only works for real-valued, symmetric positive "
269 "definite matrices. However, you instantiated this class for ScalarType"
270 " = " << Teuchos::TypeNameTraits<ScalarType>::name () << ", which is a "
271 "complex-valued type. While this may be algorithmically correct if all "
272 "of the complex numbers in the matrix have zero imaginary part, we "
273 "forbid using complex ScalarType altogether in order to remind you of "
274 "the limitations of our implementation (and of the algorithm itself).");
275 }
276 else {
277 checkInputMatrix ();
278 }
279}
280
281template<class ScalarType, class MV>
283Chebyshev (Teuchos::RCP<const row_matrix_type> A) :
284 A_ (A),
285 savedDiagOffsets_ (false),
286 computedLambdaMax_ (STS::nan ()),
287 computedLambdaMin_ (STS::nan ()),
288 lambdaMaxForApply_ (STS::nan ()),
289 lambdaMinForApply_ (STS::nan ()),
290 eigRatioForApply_ (STS::nan ()),
291 userLambdaMax_ (STS::nan ()),
292 userLambdaMin_ (STS::nan ()),
293 userEigRatio_ (Teuchos::as<ST> (30)),
294 minDiagVal_ (STS::eps ()),
295 numIters_ (1),
296 eigMaxIters_ (10),
297 eigRelTolerance_(Teuchos::ScalarTraits<MT>::zero ()),
298 eigKeepVectors_(false),
299 eigenAnalysisType_("power method"),
300 eigNormalizationFreq_(1),
301 zeroStartingSolution_ (true),
302 assumeMatrixUnchanged_ (false),
303 chebyshevAlgorithm_("first"),
304 computeMaxResNorm_ (false),
305 computeSpectralRadius_(true),
306 ckUseNativeSpMV_(false),
307 debug_ (false)
309 checkConstructorInput ();
310}
311
312template<class ScalarType, class MV>
314Chebyshev (Teuchos::RCP<const row_matrix_type> A,
315 Teuchos::ParameterList& params) :
316 A_ (A),
317 savedDiagOffsets_ (false),
318 computedLambdaMax_ (STS::nan ()),
319 computedLambdaMin_ (STS::nan ()),
320 lambdaMaxForApply_ (STS::nan ()),
321 boostFactor_ (static_cast<MT> (1.1)),
322 lambdaMinForApply_ (STS::nan ()),
323 eigRatioForApply_ (STS::nan ()),
324 userLambdaMax_ (STS::nan ()),
325 userLambdaMin_ (STS::nan ()),
326 userEigRatio_ (Teuchos::as<ST> (30)),
327 minDiagVal_ (STS::eps ()),
328 numIters_ (1),
329 eigMaxIters_ (10),
330 eigRelTolerance_(Teuchos::ScalarTraits<MT>::zero ()),
331 eigKeepVectors_(false),
332 eigenAnalysisType_("power method"),
333 eigNormalizationFreq_(1),
334 zeroStartingSolution_ (true),
335 assumeMatrixUnchanged_ (false),
336 chebyshevAlgorithm_("first"),
337 computeMaxResNorm_ (false),
338 computeSpectralRadius_(true),
339 ckUseNativeSpMV_(false),
340 debug_ (false)
341{
342 checkConstructorInput ();
343 setParameters (params);
344}
345
346template<class ScalarType, class MV>
347void
349setParameters (Teuchos::ParameterList& plist)
350{
351 using Teuchos::RCP;
352 using Teuchos::rcp;
353 using Teuchos::rcp_const_cast;
354
355 // Note to developers: The logic for this method is complicated,
356 // because we want to accept Ifpack and ML parameters whenever
357 // possible, but we don't want to add their default values to the
358 // user's ParameterList. That's why we do all the isParameter()
359 // checks, instead of using the two-argument version of get()
360 // everywhere. The min and max eigenvalue parameters are also a
361 // special case, because we decide whether or not to do eigenvalue
362 // analysis based on whether the user supplied the max eigenvalue.
363
364 // Default values of all the parameters.
365 const ST defaultLambdaMax = STS::nan ();
366 const ST defaultLambdaMin = STS::nan ();
367 // 30 is Ifpack::Chebyshev's default. ML has a different default
368 // eigRatio for smoothers and the coarse-grid solve (if using
369 // Chebyshev for that). The former uses 20; the latter uses 30.
370 // We're testing smoothers here, so use 20. (However, if you give
371 // ML an Epetra matrix, it will use Ifpack for Chebyshev, in which
372 // case it would defer to Ifpack's default settings.)
373 const ST defaultEigRatio = Teuchos::as<ST> (30);
374 const MT defaultBoostFactor = static_cast<MT> (1.1);
375 const ST defaultMinDiagVal = STS::eps ();
376 const int defaultNumIters = 1;
377 const int defaultEigMaxIters = 10;
378 const MT defaultEigRelTolerance = Teuchos::ScalarTraits<MT>::zero ();
379 const bool defaultEigKeepVectors = false;
380 const int defaultEigNormalizationFreq = 1;
381 const bool defaultZeroStartingSolution = true; // Ifpack::Chebyshev default
382 const bool defaultAssumeMatrixUnchanged = false;
383 const std::string defaultChebyshevAlgorithm = "first";
384 const bool defaultComputeMaxResNorm = false;
385 const bool defaultComputeSpectralRadius = true;
386 const bool defaultCkUseNativeSpMV = false;
387 const bool defaultDebug = false;
388
389 // We'll set the instance data transactionally, after all reads
390 // from the ParameterList. That way, if any of the ParameterList
391 // reads fail (e.g., due to the wrong parameter type), we will not
392 // have left the instance data in a half-changed state.
393 RCP<const V> userInvDiagCopy; // if nonnull: deep copy of user's Vector
394 ST lambdaMax = defaultLambdaMax;
395 ST lambdaMin = defaultLambdaMin;
396 ST eigRatio = defaultEigRatio;
397 MT boostFactor = defaultBoostFactor;
398 ST minDiagVal = defaultMinDiagVal;
399 int numIters = defaultNumIters;
400 int eigMaxIters = defaultEigMaxIters;
401 MT eigRelTolerance = defaultEigRelTolerance;
402 bool eigKeepVectors = defaultEigKeepVectors;
403 int eigNormalizationFreq = defaultEigNormalizationFreq;
404 bool zeroStartingSolution = defaultZeroStartingSolution;
405 bool assumeMatrixUnchanged = defaultAssumeMatrixUnchanged;
406 std::string chebyshevAlgorithm = defaultChebyshevAlgorithm;
407 bool computeMaxResNorm = defaultComputeMaxResNorm;
408 bool computeSpectralRadius = defaultComputeSpectralRadius;
409 bool ckUseNativeSpMV = defaultCkUseNativeSpMV;
410 bool debug = defaultDebug;
411
412 // Fetch the parameters from the ParameterList. Defer all
413 // externally visible side effects until we have finished all
414 // ParameterList interaction. This makes the method satisfy the
415 // strong exception guarantee.
416
417 if (plist.isType<bool> ("debug")) {
418 debug = plist.get<bool> ("debug");
419 }
420 else if (plist.isType<int> ("debug")) {
421 const int debugInt = plist.get<bool> ("debug");
422 debug = debugInt != 0;
423 }
424
425 // Get the user-supplied inverse diagonal.
426 //
427 // Check for a raw pointer (const V* or V*), for Ifpack
428 // compatibility, as well as for RCP<const V>, RCP<V>, const V, or
429 // V. We'll make a deep copy of the vector at the end of this
430 // method anyway, so its const-ness doesn't matter. We handle the
431 // latter two cases ("const V" or "V") specially (copy them into
432 // userInvDiagCopy first, which is otherwise null at the end of the
433 // long if-then chain) to avoid an extra copy.
434
435 const char opInvDiagLabel[] = "chebyshev: operator inv diagonal";
436 if (plist.isParameter (opInvDiagLabel)) {
437 // Pointer to the user's Vector, if provided.
438 RCP<const V> userInvDiag;
439
440 if (plist.isType<const V*> (opInvDiagLabel)) {
441 const V* rawUserInvDiag =
442 plist.get<const V*> (opInvDiagLabel);
443 // Nonowning reference (we'll make a deep copy below)
444 userInvDiag = rcp (rawUserInvDiag, false);
445 }
446 else if (plist.isType<const V*> (opInvDiagLabel)) {
447 V* rawUserInvDiag = plist.get<V*> (opInvDiagLabel);
448 // Nonowning reference (we'll make a deep copy below)
449 userInvDiag = rcp (const_cast<const V*> (rawUserInvDiag), false);
450 }
451 else if (plist.isType<RCP<const V>> (opInvDiagLabel)) {
452 userInvDiag = plist.get<RCP<const V> > (opInvDiagLabel);
453 }
454 else if (plist.isType<RCP<V>> (opInvDiagLabel)) {
455 RCP<V> userInvDiagNonConst =
456 plist.get<RCP<V> > (opInvDiagLabel);
457 userInvDiag = rcp_const_cast<const V> (userInvDiagNonConst);
458 }
459 else if (plist.isType<const V> (opInvDiagLabel)) {
460 const V& userInvDiagRef = plist.get<const V> (opInvDiagLabel);
461 userInvDiagCopy = rcp (new V (userInvDiagRef, Teuchos::Copy));
462 userInvDiag = userInvDiagCopy;
463 }
464 else if (plist.isType<V> (opInvDiagLabel)) {
465 V& userInvDiagNonConstRef = plist.get<V> (opInvDiagLabel);
466 const V& userInvDiagRef = const_cast<const V&> (userInvDiagNonConstRef);
467 userInvDiagCopy = rcp (new V (userInvDiagRef, Teuchos::Copy));
468 userInvDiag = userInvDiagCopy;
469 }
470
471 // NOTE: If the user's parameter has some strange type that we
472 // didn't test above, userInvDiag might still be null. You may
473 // want to add an error test for this condition. Currently, we
474 // just say in this case that the user didn't give us a Vector.
475
476 // If we have userInvDiag but don't have a deep copy yet, make a
477 // deep copy now.
478 if (! userInvDiag.is_null () && userInvDiagCopy.is_null ()) {
479 userInvDiagCopy = rcp (new V (*userInvDiag, Teuchos::Copy));
480 }
481
482 // NOTE: userInvDiag, if provided, is a row Map version of the
483 // Vector. We don't necessarily have a range Map yet. compute()
484 // would be the proper place to compute the range Map version of
485 // userInvDiag.
486 }
487
488 // Load the kernel fuse override from the parameter list
489 if (plist.isParameter ("chebyshev: use native spmv"))
490 ckUseNativeSpMV = plist.get("chebyshev: use native spmv", ckUseNativeSpMV);
491
492 // Don't fill in defaults for the max or min eigenvalue, because
493 // this class uses the existence of those parameters to determine
494 // whether it should do eigenanalysis.
495 if (plist.isParameter ("chebyshev: max eigenvalue")) {
496 if (plist.isType<double>("chebyshev: max eigenvalue"))
497 lambdaMax = plist.get<double> ("chebyshev: max eigenvalue");
498 else
499 lambdaMax = plist.get<ST> ("chebyshev: max eigenvalue");
500 TEUCHOS_TEST_FOR_EXCEPTION(
501 STS::isnaninf (lambdaMax), std::invalid_argument,
502 "Ifpack2::Chebyshev::setParameters: \"chebyshev: max eigenvalue\" "
503 "parameter is NaN or Inf. This parameter is optional, but if you "
504 "choose to supply it, it must have a finite value.");
505 }
506 if (plist.isParameter ("chebyshev: min eigenvalue")) {
507 if (plist.isType<double>("chebyshev: min eigenvalue"))
508 lambdaMin = plist.get<double> ("chebyshev: min eigenvalue");
509 else
510 lambdaMin = plist.get<ST> ("chebyshev: min eigenvalue");
511 TEUCHOS_TEST_FOR_EXCEPTION(
512 STS::isnaninf (lambdaMin), std::invalid_argument,
513 "Ifpack2::Chebyshev::setParameters: \"chebyshev: min eigenvalue\" "
514 "parameter is NaN or Inf. This parameter is optional, but if you "
515 "choose to supply it, it must have a finite value.");
516 }
517
518 // Only fill in Ifpack2's name for the default parameter, not ML's.
519 if (plist.isParameter ("smoother: Chebyshev alpha")) { // ML compatibility
520 if (plist.isType<double>("smoother: Chebyshev alpha"))
521 eigRatio = plist.get<double> ("smoother: Chebyshev alpha");
522 else
523 eigRatio = plist.get<ST> ("smoother: Chebyshev alpha");
524 }
525 // Ifpack2's name overrides ML's name.
526 eigRatio = plist.get ("chebyshev: ratio eigenvalue", eigRatio);
527 TEUCHOS_TEST_FOR_EXCEPTION(
528 STS::isnaninf (eigRatio), std::invalid_argument,
529 "Ifpack2::Chebyshev::setParameters: \"chebyshev: ratio eigenvalue\" "
530 "parameter (also called \"smoother: Chebyshev alpha\") is NaN or Inf. "
531 "This parameter is optional, but if you choose to supply it, it must have "
532 "a finite value.");
533 // mfh 11 Feb 2013: This class is currently only correct for real
534 // Scalar types, but we still want it to build for complex Scalar
535 // type so that users of Ifpack2::Factory can build their
536 // executables for real or complex Scalar type. Thus, we take the
537 // real parts here, which are always less-than comparable.
538 TEUCHOS_TEST_FOR_EXCEPTION(
539 STS::real (eigRatio) < STS::real (STS::one ()),
540 std::invalid_argument,
541 "Ifpack2::Chebyshev::setParameters: \"chebyshev: ratio eigenvalue\""
542 "parameter (also called \"smoother: Chebyshev alpha\") must be >= 1, "
543 "but you supplied the value " << eigRatio << ".");
544
545 // See Github Issue #234. This parameter may be either MT
546 // (preferred) or double. We check both.
547 {
548 const char paramName[] = "chebyshev: boost factor";
549
550 if (plist.isParameter (paramName)) {
551 if (plist.isType<MT> (paramName)) { // MT preferred
552 boostFactor = plist.get<MT> (paramName);
553 }
554 else if (! std::is_same<double, MT>::value &&
555 plist.isType<double> (paramName)) {
556 const double dblBF = plist.get<double> (paramName);
557 boostFactor = static_cast<MT> (dblBF);
558 }
559 else {
560 TEUCHOS_TEST_FOR_EXCEPTION
561 (true, std::invalid_argument,
562 "Ifpack2::Chebyshev::setParameters: \"chebyshev: boost factor\""
563 "parameter must have type magnitude_type (MT) or double.");
564 }
565 }
566 else { // parameter not in the list
567 // mfh 12 Aug 2016: To preserve current behavior (that fills in
568 // any parameters not in the input ParameterList with their
569 // default values), we call set() here. I don't actually like
570 // this behavior; I prefer the Belos model, where the input
571 // ParameterList is a delta from current behavior. However, I
572 // don't want to break things.
573 plist.set (paramName, defaultBoostFactor);
574 }
575 TEUCHOS_TEST_FOR_EXCEPTION
576 (boostFactor < Teuchos::ScalarTraits<MT>::one (), std::invalid_argument,
577 "Ifpack2::Chebyshev::setParameters: \"" << paramName << "\" parameter "
578 "must be >= 1, but you supplied the value " << boostFactor << ".");
579 }
580
581 // Same name in Ifpack2 and Ifpack.
582 minDiagVal = plist.get ("chebyshev: min diagonal value", minDiagVal);
583 TEUCHOS_TEST_FOR_EXCEPTION(
584 STS::isnaninf (minDiagVal), std::invalid_argument,
585 "Ifpack2::Chebyshev::setParameters: \"chebyshev: min diagonal value\" "
586 "parameter is NaN or Inf. This parameter is optional, but if you choose "
587 "to supply it, it must have a finite value.");
588
589 // Only fill in Ifpack2's name, not ML's or Ifpack's.
590 if (plist.isParameter ("smoother: sweeps")) { // ML compatibility
591 numIters = plist.get<int> ("smoother: sweeps");
592 } // Ifpack's name overrides ML's name.
593 if (plist.isParameter ("relaxation: sweeps")) { // Ifpack compatibility
594 numIters = plist.get<int> ("relaxation: sweeps");
595 } // Ifpack2's name overrides Ifpack's name.
596 numIters = plist.get ("chebyshev: degree", numIters);
597 TEUCHOS_TEST_FOR_EXCEPTION(
598 numIters < 0, std::invalid_argument,
599 "Ifpack2::Chebyshev::setParameters: \"chebyshev: degree\" parameter (also "
600 "called \"smoother: sweeps\" or \"relaxation: sweeps\") must be a "
601 "nonnegative integer. You gave a value of " << numIters << ".");
602
603 // The last parameter name overrides the first.
604 if (plist.isParameter ("eigen-analysis: iterations")) { // ML compatibility
605 eigMaxIters = plist.get<int> ("eigen-analysis: iterations");
606 } // Ifpack2's name overrides ML's name.
607 eigMaxIters = plist.get ("chebyshev: eigenvalue max iterations", eigMaxIters);
608 TEUCHOS_TEST_FOR_EXCEPTION(
609 eigMaxIters < 0, std::invalid_argument,
610 "Ifpack2::Chebyshev::setParameters: \"chebyshev: eigenvalue max iterations"
611 "\" parameter (also called \"eigen-analysis: iterations\") must be a "
612 "nonnegative integer. You gave a value of " << eigMaxIters << ".");
613
614 if (plist.isType<double>("chebyshev: eigenvalue relative tolerance"))
615 eigRelTolerance = Teuchos::as<MT>(plist.get<double> ("chebyshev: eigenvalue relative tolerance"));
616 else if (plist.isType<MT>("chebyshev: eigenvalue relative tolerance"))
617 eigRelTolerance = plist.get<MT> ("chebyshev: eigenvalue relative tolerance");
618 else if (plist.isType<ST>("chebyshev: eigenvalue relative tolerance"))
619 eigRelTolerance = Teuchos::ScalarTraits<ST>::magnitude(plist.get<ST> ("chebyshev: eigenvalue relative tolerance"));
620
621 eigKeepVectors = plist.get ("chebyshev: eigenvalue keep vectors", eigKeepVectors);
622
623 eigNormalizationFreq = plist.get ("chebyshev: eigenvalue normalization frequency", eigNormalizationFreq);
624 TEUCHOS_TEST_FOR_EXCEPTION(
625 eigNormalizationFreq < 0, std::invalid_argument,
626 "Ifpack2::Chebyshev::setParameters: \"chebyshev: eigenvalue normalization frequency"
627 "\" parameter must be a "
628 "nonnegative integer. You gave a value of " << eigNormalizationFreq << ".")
629
630 zeroStartingSolution = plist.get ("chebyshev: zero starting solution",
631 zeroStartingSolution);
632 assumeMatrixUnchanged = plist.get ("chebyshev: assume matrix does not change",
633 assumeMatrixUnchanged);
634
635 // We don't want to fill these parameters in, because they shouldn't
636 // be visible to Ifpack2::Chebyshev users.
637 if (plist.isParameter ("chebyshev: algorithm")) {
638 chebyshevAlgorithm = plist.get<std::string> ("chebyshev: algorithm");
639 TEUCHOS_TEST_FOR_EXCEPTION(
640 chebyshevAlgorithm != "first" &&
641 chebyshevAlgorithm != "textbook" &&
642 chebyshevAlgorithm != "fourth" &&
643 chebyshevAlgorithm != "opt_fourth",
644 std::invalid_argument,
645 "Ifpack2::Chebyshev: Ifpack2 only supports \"first\", \"textbook\", \"fourth\", and \"opt_fourth\", for \"chebyshev: algorithm\".");
646 }
647
648#ifdef IFPACK2_ENABLE_DEPRECATED_CODE
649 // to preserve behavior with previous input decks, only read "chebyshev:textbook algorithm" setting
650 // if a user has not specified "chebyshev: algorithm"
651 if (!plist.isParameter ("chebyshev: algorithm")) {
652 if (plist.isParameter ("chebyshev: textbook algorithm")) {
653 const bool textbookAlgorithm = plist.get<bool> ("chebyshev: textbook algorithm");
654 if(textbookAlgorithm){
655 chebyshevAlgorithm = "textbook";
656 } else {
657 chebyshevAlgorithm = "first";
658 }
659 }
660 }
661#endif
662
663 if (plist.isParameter ("chebyshev: compute max residual norm")) {
664 computeMaxResNorm = plist.get<bool> ("chebyshev: compute max residual norm");
665 }
666 if (plist.isParameter ("chebyshev: compute spectral radius")) {
667 computeSpectralRadius = plist.get<bool> ("chebyshev: compute spectral radius");
668 }
669
670
671
672 // Test for Ifpack parameters that we won't ever implement here.
673 // Be careful to use the one-argument version of get(), since the
674 // two-argment version adds the parameter if it's not there.
675 TEUCHOS_TEST_FOR_EXCEPTION
676 (plist.isType<bool> ("chebyshev: use block mode") &&
677 ! plist.get<bool> ("chebyshev: use block mode"),
678 std::invalid_argument,
679 "Ifpack2::Chebyshev requires that if you set \"chebyshev: use "
680 "block mode\" at all, you must set it to false. "
681 "Ifpack2::Chebyshev does not implement Ifpack's block mode.");
682 TEUCHOS_TEST_FOR_EXCEPTION
683 (plist.isType<bool> ("chebyshev: solve normal equations") &&
684 ! plist.get<bool> ("chebyshev: solve normal equations"),
685 std::invalid_argument,
686 "Ifpack2::Chebyshev does not and will never implement the Ifpack "
687 "parameter \"chebyshev: solve normal equations\". If you want to "
688 "solve the normal equations, construct a Tpetra::Operator that "
689 "implements A^* A, and use Chebyshev to solve A^* A x = A^* b.");
690
691 // Test for Ifpack parameters that we haven't implemented yet.
692 //
693 // For now, we only check that this ML parameter, if provided, has
694 // the one value that we support. We consider other values "invalid
695 // arguments" rather than "logic errors," because Ifpack does not
696 // implement eigenanalyses other than the power method.
697 std::string eigenAnalysisType ("power-method");
698 if (plist.isParameter ("eigen-analysis: type")) {
699 eigenAnalysisType = plist.get<std::string> ("eigen-analysis: type");
700 TEUCHOS_TEST_FOR_EXCEPTION(
701 eigenAnalysisType != "power-method" &&
702 eigenAnalysisType != "power method" &&
703 eigenAnalysisType != "cg",
704 std::invalid_argument,
705 "Ifpack2::Chebyshev: Ifpack2 only supports \"power method\" and \"cg\" for \"eigen-analysis: type\".");
706 }
707
708 // We've validated all the parameters, so it's safe now to "commit" them.
709 userInvDiag_ = userInvDiagCopy;
710 userLambdaMax_ = lambdaMax;
711 userLambdaMin_ = lambdaMin;
712 userEigRatio_ = eigRatio;
713 boostFactor_ = static_cast<MT> (boostFactor);
714 minDiagVal_ = minDiagVal;
715 numIters_ = numIters;
716 eigMaxIters_ = eigMaxIters;
717 eigRelTolerance_ = eigRelTolerance;
718 eigKeepVectors_ = eigKeepVectors;
719 eigNormalizationFreq_ = eigNormalizationFreq;
720 eigenAnalysisType_ = eigenAnalysisType;
721 zeroStartingSolution_ = zeroStartingSolution;
722 assumeMatrixUnchanged_ = assumeMatrixUnchanged;
723 chebyshevAlgorithm_ = chebyshevAlgorithm;
724 computeMaxResNorm_ = computeMaxResNorm;
725 computeSpectralRadius_ = computeSpectralRadius;
726 ckUseNativeSpMV_ = ckUseNativeSpMV;
727 debug_ = debug;
728
729 if (debug_) {
730 // Only print if myRank == 0.
731 int myRank = -1;
732 if (A_.is_null () || A_->getComm ().is_null ()) {
733 // We don't have a communicator (yet), so we assume that
734 // everybody can print. Revise this expectation in setMatrix().
735 myRank = 0;
736 }
737 else {
738 myRank = A_->getComm ()->getRank ();
739 }
740
741 if (myRank == 0) {
742 out_ = Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
743 }
744 else {
745 using Teuchos::oblackholestream; // prints nothing
746 RCP<oblackholestream> blackHole (new oblackholestream ());
747 out_ = Teuchos::getFancyOStream (blackHole);
748 }
749 }
750 else { // NOT debug
751 // free the "old" output stream, if there was one
752 out_ = Teuchos::null;
753 }
754}
755
756
757template<class ScalarType, class MV>
758void
759Chebyshev<ScalarType, MV>::reset ()
760{
761 ck_ = Teuchos::null;
762 D_ = Teuchos::null;
763 diagOffsets_ = offsets_type ();
764 savedDiagOffsets_ = false;
765 W_ = Teuchos::null;
766 computedLambdaMax_ = STS::nan ();
767 computedLambdaMin_ = STS::nan ();
768 eigVector_ = Teuchos::null;
769 eigVector2_ = Teuchos::null;
770}
771
772
773template<class ScalarType, class MV>
774void
776setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
777{
778 if (A.getRawPtr () != A_.getRawPtr ()) {
779 if (! assumeMatrixUnchanged_) {
780 reset ();
781 }
782 A_ = A;
783 ck_ = Teuchos::null; // constructed on demand
784
785 // The communicator may have changed, or we may not have had a
786 // communicator before. Thus, we may have to reset the debug
787 // output stream.
788 if (debug_) {
789 // Only print if myRank == 0.
790 int myRank = -1;
791 if (A.is_null () || A->getComm ().is_null ()) {
792 // We don't have a communicator (yet), so we assume that
793 // everybody can print. Revise this expectation in setMatrix().
794 myRank = 0;
795 }
796 else {
797 myRank = A->getComm ()->getRank ();
798 }
799
800 if (myRank == 0) {
801 out_ = Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cerr));
802 }
803 else {
804 Teuchos::RCP<Teuchos::oblackholestream> blackHole (new Teuchos::oblackholestream ());
805 out_ = Teuchos::getFancyOStream (blackHole); // print nothing on other processes
806 }
807 }
808 else { // NOT debug
809 // free the "old" output stream, if there was one
810 out_ = Teuchos::null;
811 }
812 }
813}
814
815template<class ScalarType, class MV>
816void
818{
819 using std::endl;
820 // Some of the optimizations below only work if A_ is a
821 // Tpetra::CrsMatrix. We'll make our best guess about its type
822 // here, since we have no way to get back the original fifth
823 // template parameter.
824 typedef Tpetra::CrsMatrix<typename MV::scalar_type,
825 typename MV::local_ordinal_type,
826 typename MV::global_ordinal_type,
827 typename MV::node_type> crs_matrix_type;
828
829 TEUCHOS_TEST_FOR_EXCEPTION(
830 A_.is_null (), std::runtime_error, "Ifpack2::Chebyshev::compute: The input "
831 "matrix A is null. Please call setMatrix() with a nonnull input matrix "
832 "before calling this method.");
833
834 // If A_ is a CrsMatrix and its graph is constant, we presume that
835 // the user plans to reuse the structure of A_, but possibly change
836 // A_'s values before each compute() call. This is the intended use
837 // case for caching the offsets of the diagonal entries of A_, to
838 // speed up extraction of diagonal entries on subsequent compute()
839 // calls.
840
841 // FIXME (mfh 22 Jan 2013, 10 Feb 2013) In all cases when we use
842 // isnaninf() in this method, we really only want to check if the
843 // number is NaN. Inf means something different. However,
844 // Teuchos::ScalarTraits doesn't distinguish the two cases.
845
846 // makeInverseDiagonal() returns a range Map Vector.
847 if (userInvDiag_.is_null ()) {
848 Teuchos::RCP<const crs_matrix_type> A_crsMat =
849 Teuchos::rcp_dynamic_cast<const crs_matrix_type> (A_);
850 if (D_.is_null ()) { // We haven't computed D_ before
851 if (! A_crsMat.is_null () && A_crsMat->isFillComplete ()) {
852 // It's a CrsMatrix with a const graph; cache diagonal offsets.
853 const size_t lclNumRows = A_crsMat->getLocalNumRows ();
854 if (diagOffsets_.extent (0) < lclNumRows) {
855 diagOffsets_ = offsets_type (); // clear first to save memory
856 diagOffsets_ = offsets_type ("offsets", lclNumRows);
857 }
858 A_crsMat->getCrsGraph ()->getLocalDiagOffsets (diagOffsets_);
859 savedDiagOffsets_ = true;
860 D_ = makeInverseDiagonal (*A_, true);
861 }
862 else { // either A_ is not a CrsMatrix, or its graph is nonconst
863 D_ = makeInverseDiagonal (*A_);
864 }
865 }
866 else if (! assumeMatrixUnchanged_) { // D_ exists but A_ may have changed
867 if (! A_crsMat.is_null () && A_crsMat->isFillComplete ()) {
868 // It's a CrsMatrix with a const graph; cache diagonal offsets
869 // if we haven't already.
870 if (! savedDiagOffsets_) {
871 const size_t lclNumRows = A_crsMat->getLocalNumRows ();
872 if (diagOffsets_.extent (0) < lclNumRows) {
873 diagOffsets_ = offsets_type (); // clear first to save memory
874 diagOffsets_ = offsets_type ("offsets", lclNumRows);
875 }
876 A_crsMat->getCrsGraph ()->getLocalDiagOffsets (diagOffsets_);
877 savedDiagOffsets_ = true;
878 }
879 // Now we're guaranteed to have cached diagonal offsets.
880 D_ = makeInverseDiagonal (*A_, true);
881 }
882 else { // either A_ is not a CrsMatrix, or its graph is nonconst
883 D_ = makeInverseDiagonal (*A_);
884 }
885 }
886 }
887 else { // the user provided an inverse diagonal
888 D_ = makeRangeMapVectorConst (userInvDiag_);
889 }
890
891 // Have we estimated eigenvalues before?
892 const bool computedEigenvalueEstimates =
893 STS::isnaninf (computedLambdaMax_) || STS::isnaninf (computedLambdaMin_);
894
895 // Only recompute the eigenvalue estimates if
896 // - we are supposed to assume that the matrix may have changed, or
897 // - they haven't been computed before, and the user hasn't given
898 // us at least an estimate of the max eigenvalue.
899 //
900 // We at least need an estimate of the max eigenvalue. This is the
901 // most important one if using Chebyshev as a smoother.
902
903 if (! assumeMatrixUnchanged_ ||
904 (! computedEigenvalueEstimates && STS::isnaninf (userLambdaMax_))) {
905 ST computedLambdaMax;
906 if ((eigenAnalysisType_ == "power method") || (eigenAnalysisType_ == "power-method")) {
907 Teuchos::RCP<V> x;
908 if (eigVector_.is_null()) {
909 x = Teuchos::rcp(new V(A_->getDomainMap ()));
910 if (eigKeepVectors_)
911 eigVector_ = x;
913 } else
914 x = eigVector_;
915
916 Teuchos::RCP<V> y;
917 if (eigVector2_.is_null()) {
918 y = rcp(new V(A_->getRangeMap ()));
919 if (eigKeepVectors_)
920 eigVector2_ = y;
921 } else
922 y = eigVector2_;
923
924 Teuchos::RCP<Teuchos::FancyOStream> stream = (debug_ ? out_ : Teuchos::null);
925 computedLambdaMax = PowerMethod::powerMethodWithInitGuess (*A_, *D_, eigMaxIters_, x, y,
926 eigRelTolerance_, eigNormalizationFreq_, stream,
927 computeSpectralRadius_);
928 }
929 else {
930 computedLambdaMax = cgMethod (*A_, *D_, eigMaxIters_);
931 }
932 TEUCHOS_TEST_FOR_EXCEPTION(
933 STS::isnaninf (computedLambdaMax),
934 std::runtime_error,
935 "Ifpack2::Chebyshev::compute: Estimation of the max eigenvalue "
936 "of D^{-1} A failed, by producing Inf or NaN. This probably means that "
937 "the matrix contains Inf or NaN values, or that it is badly scaled.");
938 TEUCHOS_TEST_FOR_EXCEPTION(
939 STS::isnaninf (userEigRatio_),
940 std::logic_error,
941 "Ifpack2::Chebyshev::compute: userEigRatio_ is Inf or NaN."
942 << endl << "This should be impossible." << endl <<
943 "Please report this bug to the Ifpack2 developers.");
944
945 // The power method doesn't estimate the min eigenvalue, so we
946 // do our best to provide an estimate. userEigRatio_ has a
947 // reasonable default value, and if the user provided it, we
948 // have already checked that its value is finite and >= 1.
949 const ST computedLambdaMin = computedLambdaMax / userEigRatio_;
950
951 // Defer "committing" results until all computations succeeded.
952 computedLambdaMax_ = computedLambdaMax;
953 computedLambdaMin_ = computedLambdaMin;
954 } else {
955 TEUCHOS_TEST_FOR_EXCEPTION(
956 STS::isnaninf (userLambdaMax_) && STS::isnaninf (computedLambdaMax_),
957 std::logic_error,
958 "Ifpack2::Chebyshev::compute: " << endl <<
959 "Both userLambdaMax_ and computedLambdaMax_ are Inf or NaN."
960 << endl << "This should be impossible." << endl <<
961 "Please report this bug to the Ifpack2 developers.");
962 }
963
965 // Figure out the eigenvalue estimates that apply() will use.
967
968 // Always favor the user's max eigenvalue estimate, if provided.
969 lambdaMaxForApply_ = STS::isnaninf (userLambdaMax_) ? computedLambdaMax_ : userLambdaMax_;
970
971 // mfh 11 Feb 2013: For now, we imitate Ifpack by ignoring the
972 // user's min eigenvalue estimate, and using the given eigenvalue
973 // ratio to estimate the min eigenvalue. We could instead do this:
974 // favor the user's eigenvalue ratio estimate, but if it's not
975 // provided, use lambdaMax / lambdaMin. However, we want Chebyshev
976 // to have sensible smoother behavior if the user did not provide
977 // eigenvalue estimates. Ifpack's behavior attempts to push down
978 // the error terms associated with the largest eigenvalues, while
979 // expecting that users will only want a small number of iterations,
980 // so that error terms associated with the smallest eigenvalues
981 // won't grow too much. This is sensible behavior for a smoother.
982 lambdaMinForApply_ = lambdaMaxForApply_ / userEigRatio_;
983 eigRatioForApply_ = userEigRatio_;
984
985 if (chebyshevAlgorithm_ == "first") {
986 // Ifpack has a special-case modification of the eigenvalue bounds
987 // for the case where the max eigenvalue estimate is close to one.
988 const ST one = Teuchos::as<ST> (1);
989 // FIXME (mfh 20 Nov 2013) Should scale this 1.0e-6 term
990 // appropriately for MT's machine precision.
991 if (STS::magnitude (lambdaMaxForApply_ - one) < Teuchos::as<MT> (1.0e-6)) {
992 lambdaMinForApply_ = one;
993 lambdaMaxForApply_ = lambdaMinForApply_;
994 eigRatioForApply_ = one; // Ifpack doesn't include this line.
995 }
996 }
997}
998
999
1000template<class ScalarType, class MV>
1001ScalarType
1003getLambdaMaxForApply() const {
1004 return lambdaMaxForApply_;
1005}
1006
1007
1008template<class ScalarType, class MV>
1011{
1012 const char prefix[] = "Ifpack2::Chebyshev::apply: ";
1013
1014 if (debug_) {
1015 *out_ << "apply: " << std::endl;
1016 }
1017 TEUCHOS_TEST_FOR_EXCEPTION
1018 (A_.is_null (), std::runtime_error, prefix << "The input matrix A is null. "
1019 " Please call setMatrix() with a nonnull input matrix, and then call "
1020 "compute(), before calling this method.");
1021 TEUCHOS_TEST_FOR_EXCEPTION
1022 (STS::isnaninf (lambdaMaxForApply_), std::runtime_error,
1023 prefix << "There is no estimate for the max eigenvalue."
1024 << std::endl << computeBeforeApplyReminder);
1025 TEUCHOS_TEST_FOR_EXCEPTION
1026 (STS::isnaninf (lambdaMinForApply_), std::runtime_error,
1027 prefix << "There is no estimate for the min eigenvalue."
1028 << std::endl << computeBeforeApplyReminder);
1029 TEUCHOS_TEST_FOR_EXCEPTION
1030 (STS::isnaninf (eigRatioForApply_), std::runtime_error,
1031 prefix <<"There is no estimate for the ratio of the max "
1032 "eigenvalue to the min eigenvalue."
1033 << std::endl << computeBeforeApplyReminder);
1034 TEUCHOS_TEST_FOR_EXCEPTION
1035 (D_.is_null (), std::runtime_error, prefix << "The vector of inverse "
1036 "diagonal entries of the matrix has not yet been computed."
1037 << std::endl << computeBeforeApplyReminder);
1038
1039 if (chebyshevAlgorithm_ == "fourth" || chebyshevAlgorithm_ == "opt_fourth") {
1040 fourthKindApplyImpl (*A_, B, X, numIters_, lambdaMaxForApply_, *D_);
1041 }
1042 else if (chebyshevAlgorithm_ == "textbook") {
1043 textbookApplyImpl (*A_, B, X, numIters_, lambdaMaxForApply_,
1044 lambdaMinForApply_, eigRatioForApply_, *D_);
1045 }
1046 else {
1047 ifpackApplyImpl (*A_, B, X, numIters_, lambdaMaxForApply_,
1048 lambdaMinForApply_, eigRatioForApply_, *D_);
1049 }
1050
1051 if (computeMaxResNorm_ && B.getNumVectors () > 0) {
1052 MV R (B.getMap (), B.getNumVectors ());
1053 computeResidual (R, B, *A_, X);
1054 Teuchos::Array<MT> norms (B.getNumVectors ());
1055 R.norm2 (norms ());
1056 return *std::max_element (norms.begin (), norms.end ());
1057 }
1058 else {
1059 return Teuchos::ScalarTraits<MT>::zero ();
1060 }
1061}
1062
1063template<class ScalarType, class MV>
1064void
1066print (std::ostream& out)
1067{
1068 using Teuchos::rcpFromRef;
1069 this->describe (* (Teuchos::getFancyOStream (rcpFromRef (out))),
1070 Teuchos::VERB_MEDIUM);
1071}
1072
1073template<class ScalarType, class MV>
1074void
1077 const ScalarType& alpha,
1078 const V& D_inv,
1079 const MV& B,
1080 MV& X)
1081{
1082 solve (W, alpha, D_inv, B); // W = alpha*D_inv*B
1083 Tpetra::deep_copy (X, W); // X = 0 + W
1084}
1085
1086template<class ScalarType, class MV>
1087void
1089computeResidual (MV& R, const MV& B, const op_type& A, const MV& X,
1090 const Teuchos::ETransp mode)
1091{
1092 Tpetra::Details::residual(A,X,B,R);
1093}
1094
1095template <class ScalarType, class MV>
1096void
1097Chebyshev<ScalarType, MV>::
1098solve (MV& Z, const V& D_inv, const MV& R) {
1099 Z.elementWiseMultiply (STS::one(), D_inv, R, STS::zero());
1100}
1101
1102template<class ScalarType, class MV>
1103void
1105solve (MV& Z, const ST alpha, const V& D_inv, const MV& R) {
1106 Z.elementWiseMultiply (alpha, D_inv, R, STS::zero());
1107}
1108
1109template<class ScalarType, class MV>
1110Teuchos::RCP<const typename Chebyshev<ScalarType, MV>::V>
1112makeInverseDiagonal (const row_matrix_type& A, const bool useDiagOffsets) const
1113{
1114 using Teuchos::RCP;
1115 using Teuchos::rcpFromRef;
1116 using Teuchos::rcp_dynamic_cast;
1117
1118 RCP<V> D_rowMap;
1119 if (!D_.is_null() &&
1120 D_->getMap()->isSameAs(*(A.getRowMap ()))) {
1121 if (debug_)
1122 *out_ << "Reusing pre-existing vector for diagonal extraction" << std::endl;
1123 D_rowMap = Teuchos::rcp_const_cast<V>(D_);
1124 } else {
1125 D_rowMap = Teuchos::rcp(new V (A.getRowMap (), /*zeroOut=*/false));
1126 if (debug_)
1127 *out_ << "Allocated new vector for diagonal extraction" << std::endl;
1128 }
1129 if (useDiagOffsets) {
1130 // The optimizations below only work if A_ is a Tpetra::CrsMatrix.
1131 // We'll make our best guess about its type here, since we have no
1132 // way to get back the original fifth template parameter.
1133 typedef Tpetra::CrsMatrix<typename MV::scalar_type,
1134 typename MV::local_ordinal_type,
1135 typename MV::global_ordinal_type,
1136 typename MV::node_type> crs_matrix_type;
1137 RCP<const crs_matrix_type> A_crsMat =
1138 rcp_dynamic_cast<const crs_matrix_type> (rcpFromRef (A));
1139 if (! A_crsMat.is_null ()) {
1140 TEUCHOS_TEST_FOR_EXCEPTION(
1141 ! savedDiagOffsets_, std::logic_error,
1142 "Ifpack2::Details::Chebyshev::makeInverseDiagonal: "
1143 "It is not allowed to call this method with useDiagOffsets=true, "
1144 "if you have not previously saved offsets of diagonal entries. "
1145 "This situation should never arise if this class is used properly. "
1146 "Please report this bug to the Ifpack2 developers.");
1147 A_crsMat->getLocalDiagCopy (*D_rowMap, diagOffsets_);
1148 }
1149 }
1150 else {
1151 // This always works for a Tpetra::RowMatrix, even if it is not a
1152 // Tpetra::CrsMatrix. We just don't have offsets in this case.
1153 A.getLocalDiagCopy (*D_rowMap);
1154 }
1155 RCP<V> D_rangeMap = makeRangeMapVector (D_rowMap);
1156
1157 if (debug_) {
1158 // In debug mode, make sure that all diagonal entries are
1159 // positive, on all processes. Note that *out_ only prints on
1160 // Process 0 of the matrix's communicator.
1161 bool foundNonpositiveValue = false;
1162 {
1163 auto D_lcl = D_rangeMap->getLocalViewHost (Tpetra::Access::ReadOnly);
1164 auto D_lcl_1d = Kokkos::subview (D_lcl, Kokkos::ALL (), 0);
1165
1166 typedef typename MV::impl_scalar_type IST;
1167 typedef typename MV::local_ordinal_type LO;
1168 typedef Kokkos::ArithTraits<IST> ATS;
1169 typedef Kokkos::ArithTraits<typename ATS::mag_type> STM;
1170
1171 const LO lclNumRows = static_cast<LO> (D_rangeMap->getLocalLength ());
1172 for (LO i = 0; i < lclNumRows; ++i) {
1173 if (STS::real (D_lcl_1d(i)) <= STM::zero ()) {
1174 foundNonpositiveValue = true;
1175 break;
1176 }
1177 }
1178 }
1179
1180 using Teuchos::outArg;
1181 using Teuchos::REDUCE_MIN;
1182 using Teuchos::reduceAll;
1183
1184 const int lclSuccess = foundNonpositiveValue ? 0 : 1;
1185 int gblSuccess = lclSuccess; // to be overwritten
1186 if (! D_rangeMap->getMap ().is_null () && D_rangeMap->getMap ()->getComm ().is_null ()) {
1187 const Teuchos::Comm<int>& comm = * (D_rangeMap->getMap ()->getComm ());
1188 reduceAll<int, int> (comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
1189 }
1190 if (gblSuccess == 1) {
1191 *out_ << "makeInverseDiagonal: The matrix's diagonal entries all have "
1192 "positive real part (this is good for Chebyshev)." << std::endl;
1193 }
1194 else {
1195 *out_ << "makeInverseDiagonal: The matrix's diagonal has at least one "
1196 "entry with nonpositive real part, on at least one process in the "
1197 "matrix's communicator. This is bad for Chebyshev." << std::endl;
1198 }
1199 }
1200
1201 // Invert the diagonal entries, replacing entries less (in
1202 // magnitude) than the user-specified value with that value.
1203 reciprocal_threshold (*D_rangeMap, minDiagVal_);
1204 return Teuchos::rcp_const_cast<const V> (D_rangeMap);
1205}
1206
1207
1208template<class ScalarType, class MV>
1209Teuchos::RCP<const typename Chebyshev<ScalarType, MV>::V>
1211makeRangeMapVectorConst (const Teuchos::RCP<const V>& D) const
1212{
1213 using Teuchos::RCP;
1214 using Teuchos::rcp;
1215 typedef Tpetra::Export<typename MV::local_ordinal_type,
1216 typename MV::global_ordinal_type,
1217 typename MV::node_type> export_type;
1218 // This throws logic_error instead of runtime_error, because the
1219 // methods that call makeRangeMapVector should all have checked
1220 // whether A_ is null before calling this method.
1221 TEUCHOS_TEST_FOR_EXCEPTION(
1222 A_.is_null (), std::logic_error, "Ifpack2::Details::Chebyshev::"
1223 "makeRangeMapVector: The input matrix A is null. Please call setMatrix() "
1224 "with a nonnull input matrix before calling this method. This is probably "
1225 "a bug in Ifpack2; please report this bug to the Ifpack2 developers.");
1226 TEUCHOS_TEST_FOR_EXCEPTION(
1227 D.is_null (), std::logic_error, "Ifpack2::Details::Chebyshev::"
1228 "makeRangeMapVector: The input Vector D is null. This is probably "
1229 "a bug in Ifpack2; please report this bug to the Ifpack2 developers.");
1230
1231 RCP<const map_type> sourceMap = D->getMap ();
1232 RCP<const map_type> rangeMap = A_->getRangeMap ();
1233 RCP<const map_type> rowMap = A_->getRowMap ();
1234
1235 if (rangeMap->isSameAs (*sourceMap)) {
1236 // The given vector's Map is the same as the matrix's range Map.
1237 // That means we don't need to Export. This is the fast path.
1238 return D;
1239 }
1240 else { // We need to Export.
1241 RCP<const export_type> exporter;
1242 // Making an Export object from scratch is expensive enough that
1243 // it's worth the O(1) global reductions to call isSameAs(), to
1244 // see if we can avoid that cost.
1245 if (sourceMap->isSameAs (*rowMap)) {
1246 // We can reuse the matrix's Export object, if there is one.
1247 exporter = A_->getGraph ()->getExporter ();
1248 }
1249 else { // We have to make a new Export object.
1250 exporter = rcp (new export_type (sourceMap, rangeMap));
1251 }
1252
1253 if (exporter.is_null ()) {
1254 return D; // Row Map and range Map are the same; no need to Export.
1255 }
1256 else { // Row Map and range Map are _not_ the same; must Export.
1257 RCP<V> D_out = rcp (new V (*D, Teuchos::Copy));
1258 D_out->doExport (*D, *exporter, Tpetra::ADD);
1259 return Teuchos::rcp_const_cast<const V> (D_out);
1260 }
1261 }
1262}
1263
1264
1265template<class ScalarType, class MV>
1266Teuchos::RCP<typename Chebyshev<ScalarType, MV>::V>
1268makeRangeMapVector (const Teuchos::RCP<V>& D) const
1269{
1270 using Teuchos::rcp_const_cast;
1271 return rcp_const_cast<V> (makeRangeMapVectorConst (rcp_const_cast<V> (D)));
1272}
1273
1274
1275template<class ScalarType, class MV>
1276void
1278textbookApplyImpl (const op_type& A,
1279 const MV& B,
1280 MV& X,
1281 const int numIters,
1282 const ST lambdaMax,
1283 const ST lambdaMin,
1284 const ST eigRatio,
1285 const V& D_inv) const
1286{
1287 (void) lambdaMin; // Forestall compiler warning.
1288 const ST myLambdaMin = lambdaMax / eigRatio;
1289
1290 const ST zero = Teuchos::as<ST> (0);
1291 const ST one = Teuchos::as<ST> (1);
1292 const ST two = Teuchos::as<ST> (2);
1293 const ST d = (lambdaMax + myLambdaMin) / two; // Ifpack2 calls this theta
1294 const ST c = (lambdaMax - myLambdaMin) / two; // Ifpack2 calls this 1/delta
1295
1296 if (zeroStartingSolution_ && numIters > 0) {
1297 // If zero iterations, then input X is output X.
1298 X.putScalar (zero);
1299 }
1300 MV R (B.getMap (), B.getNumVectors (), false);
1301 MV P (B.getMap (), B.getNumVectors (), false);
1302 MV Z (B.getMap (), B.getNumVectors (), false);
1303 ST alpha, beta;
1304 for (int i = 0; i < numIters; ++i) {
1305 computeResidual (R, B, A, X); // R = B - A*X
1306 solve (Z, D_inv, R); // z = D_inv * R, that is, D \ R.
1307 if (i == 0) {
1308 P = Z;
1309 alpha = two / d;
1310 } else {
1311 //beta = (c * alpha / two)^2;
1312 //const ST sqrtBeta = c * alpha / two;
1313 //beta = sqrtBeta * sqrtBeta;
1314 beta = alpha * (c/two) * (c/two);
1315 alpha = one / (d - beta);
1316 P.update (one, Z, beta); // P = Z + beta*P
1317 }
1318 X.update (alpha, P, one); // X = X + alpha*P
1319 // If we compute the residual here, we could either do R = B -
1320 // A*X, or R = R - alpha*A*P. Since we choose the former, we
1321 // can move the computeResidual call to the top of the loop.
1322 }
1323}
1324
1325template<class ScalarType, class MV>
1326void
1328fourthKindApplyImpl (const op_type& A,
1329 const MV& B,
1330 MV& X,
1331 const int numIters,
1332 const ST lambdaMax,
1333 const V& D_inv)
1334{
1335 // standard 4th kind Chebyshev smoother has \beta_i := 1
1336 std::vector<ScalarType> betas(numIters, 1.0);
1337 if(chebyshevAlgorithm_ == "opt_fourth"){
1338 betas = optimalWeightsImpl<ScalarType>(numIters);
1339 }
1340
1341 const ST invEig = MT(1) / (lambdaMax * boostFactor_);
1342
1343 // Fetch cached temporary (multi)vector.
1344 Teuchos::RCP<MV> Z_ptr = makeTempMultiVector (B);
1345 MV& Z = *Z_ptr;
1346
1347 // Store 4th-kind result (needed as temporary for bootstrapping opt. 4th-kind Chebyshev)
1348 // Fetch the second cached temporary (multi)vector.
1349 Teuchos::RCP<MV> X4_ptr = makeSecondTempMultiVector (B);
1350 MV& X4 = *X4_ptr;
1351
1352 // Special case for the first iteration.
1353 if (! zeroStartingSolution_) {
1354
1355 // X4 = X
1356 Tpetra::deep_copy (X4, X);
1357
1358 if (ck_.is_null ()) {
1359 Teuchos::RCP<const op_type> A_op = A_;
1360 ck_ = Teuchos::rcp (new ChebyshevKernel<op_type> (A_op, ckUseNativeSpMV_));
1361 }
1362 // Z := (4/3 * invEig)*D_inv*(B-A*X4)
1363 // X4 := X4 + Z
1364 ck_->compute (Z, MT(4.0/3.0) * invEig, const_cast<V&> (D_inv),
1365 const_cast<MV&> (B), X4, STS::zero());
1366
1367 // X := X + beta[0] * Z
1368 X.update (betas[0], Z, STS::one());
1369 }
1370 else {
1371 // Z := (4/3 * invEig)*D_inv*B and X := 0 + Z.
1372 firstIterationWithZeroStartingSolution (Z, MT(4.0/3.0) * invEig, D_inv, B, X4);
1373
1374 // X := 0 + beta * Z
1375 X.update (betas[0], Z, STS::zero());
1376 }
1377
1378 if (numIters > 1 && ck_.is_null ()) {
1379 Teuchos::RCP<const op_type> A_op = A_;
1380 ck_ = Teuchos::rcp (new ChebyshevKernel<op_type> (A_op, ckUseNativeSpMV_));
1381 }
1382
1383 for (int i = 1; i < numIters; ++i) {
1384 const ST zScale = (2.0 * i - 1.0) / (2.0 * i + 3.0);
1385 const ST rScale = MT((8.0 * i + 4.0) / (2.0 * i + 3.0)) * invEig;
1386
1387 // Z := rScale*D_inv*(B - A*X4) + zScale*Z.
1388 // X4 := X4 + Z
1389 ck_->compute (Z, rScale, const_cast<V&> (D_inv),
1390 const_cast<MV&> (B), (X4), zScale);
1391
1392 // X := X + beta[i] * Z
1393 X.update (betas[i], Z, STS::one());
1394 }
1395}
1396
1397template<class ScalarType, class MV>
1399Chebyshev<ScalarType, MV>::maxNormInf (const MV& X) {
1400 Teuchos::Array<MT> norms (X.getNumVectors ());
1401 X.normInf (norms());
1402 return *std::max_element (norms.begin (), norms.end ());
1403}
1404
1405template<class ScalarType, class MV>
1406void
1408ifpackApplyImpl (const op_type& A,
1409 const MV& B,
1410 MV& X,
1411 const int numIters,
1412 const ST lambdaMax,
1413 const ST lambdaMin,
1414 const ST eigRatio,
1415 const V& D_inv)
1416{
1417 using std::endl;
1418#ifdef HAVE_IFPACK2_DEBUG
1419 const bool debug = debug_;
1420#else
1421 const bool debug = false;
1422#endif
1423
1424 if (debug) {
1425 *out_ << " \\|B\\|_{\\infty} = " << maxNormInf (B) << endl;
1426 *out_ << " \\|X\\|_{\\infty} = " << maxNormInf (X) << endl;
1427 }
1428
1429 if (numIters <= 0) {
1430 return;
1431 }
1432 const ST zero = static_cast<ST> (0.0);
1433 const ST one = static_cast<ST> (1.0);
1434 const ST two = static_cast<ST> (2.0);
1435
1436 // Quick solve when the matrix A is the identity.
1437 if (lambdaMin == one && lambdaMax == lambdaMin) {
1438 solve (X, D_inv, B);
1439 return;
1440 }
1441
1442 // Initialize coefficients
1443 const ST alpha = lambdaMax / eigRatio;
1444 const ST beta = boostFactor_ * lambdaMax;
1445 const ST delta = two / (beta - alpha);
1446 const ST theta = (beta + alpha) / two;
1447 const ST s1 = theta * delta;
1448
1449 if (debug) {
1450 *out_ << " alpha = " << alpha << endl
1451 << " beta = " << beta << endl
1452 << " delta = " << delta << endl
1453 << " theta = " << theta << endl
1454 << " s1 = " << s1 << endl;
1455 }
1456
1457 // Fetch cached temporary (multi)vector.
1458 Teuchos::RCP<MV> W_ptr = makeTempMultiVector (B);
1459 MV& W = *W_ptr;
1460
1461 if (debug) {
1462 *out_ << " Iteration " << 1 << ":" << endl
1463 << " - \\|D\\|_{\\infty} = " << D_->normInf () << endl;
1464 }
1465
1466 // Special case for the first iteration.
1467 if (! zeroStartingSolution_) {
1468 // mfh 22 May 2019: Tests don't actually exercise this path.
1469
1470 if (ck_.is_null ()) {
1471 Teuchos::RCP<const op_type> A_op = A_;
1472 ck_ = Teuchos::rcp (new ChebyshevKernel<op_type> (A_op, ckUseNativeSpMV_));
1473 }
1474 // W := (1/theta)*D_inv*(B-A*X) and X := X + W.
1475 // X := X + W
1476 ck_->compute (W, one/theta, const_cast<V&> (D_inv),
1477 const_cast<MV&> (B), X, zero);
1478 }
1479 else {
1480 // W := (1/theta)*D_inv*B and X := 0 + W.
1481 firstIterationWithZeroStartingSolution (W, one/theta, D_inv, B, X);
1482 }
1483
1484 if (debug) {
1485 *out_ << " - \\|W\\|_{\\infty} = " << maxNormInf (W) << endl
1486 << " - \\|X\\|_{\\infty} = " << maxNormInf (X) << endl;
1487 }
1488
1489 if (numIters > 1 && ck_.is_null ()) {
1490 Teuchos::RCP<const op_type> A_op = A_;
1491 ck_ = Teuchos::rcp (new ChebyshevKernel<op_type> (A_op, ckUseNativeSpMV_));
1492 }
1493
1494 // The rest of the iterations.
1495 ST rhok = one / s1;
1496 ST rhokp1, dtemp1, dtemp2;
1497 for (int deg = 1; deg < numIters; ++deg) {
1498 if (debug) {
1499 *out_ << " Iteration " << deg+1 << ":" << endl
1500 << " - \\|D\\|_{\\infty} = " << D_->normInf () << endl
1501 << " - \\|B\\|_{\\infty} = " << maxNormInf (B) << endl
1502 << " - \\|A\\|_{\\text{frob}} = " << A_->getFrobeniusNorm ()
1503 << endl << " - rhok = " << rhok << endl;
1504 }
1505
1506 rhokp1 = one / (two * s1 - rhok);
1507 dtemp1 = rhokp1 * rhok;
1508 dtemp2 = two * rhokp1 * delta;
1509 rhok = rhokp1;
1510
1511 if (debug) {
1512 *out_ << " - dtemp1 = " << dtemp1 << endl
1513 << " - dtemp2 = " << dtemp2 << endl;
1514 }
1515
1516 // W := dtemp2*D_inv*(B - A*X) + dtemp1*W.
1517 // X := X + W
1518 ck_->compute (W, dtemp2, const_cast<V&> (D_inv),
1519 const_cast<MV&> (B), (X), dtemp1);
1520
1521 if (debug) {
1522 *out_ << " - \\|W\\|_{\\infty} = " << maxNormInf (W) << endl
1523 << " - \\|X\\|_{\\infty} = " << maxNormInf (X) << endl;
1524 }
1525 }
1526}
1527
1528
1529
1530
1531template<class ScalarType, class MV>
1534cgMethodWithInitGuess (const op_type& A,
1535 const V& D_inv,
1536 const int numIters,
1537 V& r)
1538{
1539 using std::endl;
1540 using MagnitudeType = typename STS::magnitudeType;
1541 if (debug_) {
1542 *out_ << " cgMethodWithInitGuess:" << endl;
1543 }
1544
1545
1546
1547 const ST one = STS::one();
1548 ST beta, betaOld = one, pAp, pApOld = one, alpha, rHz, rHzOld, rHzOld2 = one, lambdaMax;
1549 // ST lambdaMin;
1550 Teuchos::ArrayRCP<MagnitudeType> diag, offdiag;
1551 Teuchos::RCP<V> p, z, Ap;
1552 diag.resize(numIters);
1553 offdiag.resize(numIters-1);
1554
1555 p = rcp(new V(A.getRangeMap ()));
1556 z = rcp(new V(A.getRangeMap ()));
1557 Ap = rcp(new V(A.getRangeMap ()));
1558
1559 // Tpetra::Details::residual (A, x, *b, *r);
1560 solve (*p, D_inv, r);
1561 rHz = r.dot (*p);
1562
1563 for (int iter = 0; iter < numIters; ++iter) {
1564 if (debug_) {
1565 *out_ << " Iteration " << iter << endl;
1566 }
1567 A.apply (*p, *Ap);
1568 pAp = p->dot (*Ap);
1569 alpha = rHz/pAp;
1570 r.update (-alpha, *Ap, one);
1571 rHzOld = rHz;
1572 solve (*z, D_inv, r);
1573 rHz = r.dot (*z);
1574 beta = rHz / rHzOld;
1575 p->update (one, *z, beta);
1576 if (iter>0) {
1577 diag[iter] = STS::real((betaOld*betaOld * pApOld + pAp) / rHzOld);
1578 offdiag[iter-1] = -STS::real((betaOld * pApOld / (sqrt(rHzOld * rHzOld2))));
1579 if (debug_) {
1580 *out_ << " diag[" << iter << "] = " << diag[iter] << endl;
1581 *out_ << " offdiag["<< iter-1 << "] = " << offdiag[iter-1] << endl;
1582 *out_ << " rHz = "<<rHz <<endl;
1583 *out_ << " alpha = "<<alpha<<endl;
1584 *out_ << " beta = "<<beta<<endl;
1585 }
1586 } else {
1587 diag[iter] = STS::real(pAp/rHzOld);
1588 if (debug_) {
1589 *out_ << " diag[" << iter << "] = " << diag[iter] << endl;
1590 *out_ << " rHz = "<<rHz <<endl;
1591 *out_ << " alpha = "<<alpha<<endl;
1592 *out_ << " beta = "<<beta<<endl;
1593 }
1594 }
1595 rHzOld2 = rHzOld;
1596 betaOld = beta;
1597 pApOld = pAp;
1598 }
1599
1600 lambdaMax = LapackHelper<ST>::tri_diag_spectral_radius(diag, offdiag);
1601
1602 return lambdaMax;
1603}
1604
1605
1606template<class ScalarType, class MV>
1609cgMethod (const op_type& A, const V& D_inv, const int numIters)
1610{
1611 using std::endl;
1612
1613 if (debug_) {
1614 *out_ << "cgMethod:" << endl;
1615 }
1616
1617 Teuchos::RCP<V> r;
1618 if (eigVector_.is_null()) {
1619 r = rcp(new V(A.getDomainMap ()));
1620 if (eigKeepVectors_)
1621 eigVector_ = r;
1622 // For CG, we need to get the BCs right and we'll use D_inv to get that
1623 Details::computeInitialGuessForCG (D_inv,*r);
1624 } else
1625 r = eigVector_;
1626
1627 ST lambdaMax = cgMethodWithInitGuess (A, D_inv, numIters, *r);
1628
1629 return lambdaMax;
1630}
1631
1632template<class ScalarType, class MV>
1633Teuchos::RCP<const typename Chebyshev<ScalarType, MV>::row_matrix_type>
1635 return A_;
1636}
1637
1638template<class ScalarType, class MV>
1639bool
1641hasTransposeApply () const {
1642 // Technically, this is true, because the matrix must be symmetric.
1643 return true;
1644}
1645
1646template<class ScalarType, class MV>
1647Teuchos::RCP<MV>
1649makeTempMultiVector (const MV& B)
1650{
1651 // ETP 02/08/17: We must check not only if the temporary vectors are
1652 // null, but also if the number of columns match, since some multi-RHS
1653 // solvers (e.g., Belos) may call apply() with different numbers of columns.
1654
1655 const size_t B_numVecs = B.getNumVectors ();
1656 if (W_.is_null () || W_->getNumVectors () != B_numVecs) {
1657 W_ = Teuchos::rcp (new MV (B.getMap (), B_numVecs, false));
1658 }
1659 return W_;
1660}
1661
1662template<class ScalarType, class MV>
1663Teuchos::RCP<MV>
1664Chebyshev<ScalarType, MV>::
1665makeSecondTempMultiVector (const MV& B)
1666{
1667 // ETP 02/08/17: We must check not only if the temporary vectors are
1668 // null, but also if the number of columns match, since some multi-RHS
1669 // solvers (e.g., Belos) may call apply() with different numbers of columns.
1670
1671 const size_t B_numVecs = B.getNumVectors ();
1672 if (W2_.is_null () || W2_->getNumVectors () != B_numVecs) {
1673 W2_ = Teuchos::rcp (new MV (B.getMap (), B_numVecs, false));
1674 }
1675 return W2_;
1676}
1677
1678
1679template<class ScalarType, class MV>
1680std::string
1682description () const {
1683 std::ostringstream oss;
1684 // YAML requires quoting the key in this case, to distinguish
1685 // key's colons from the colon that separates key from value.
1686 oss << "\"Ifpack2::Details::Chebyshev\":"
1687 << "{"
1688 << "degree: " << numIters_
1689 << ", lambdaMax: " << lambdaMaxForApply_
1690 << ", alpha: " << eigRatioForApply_
1691 << ", lambdaMin: " << lambdaMinForApply_
1692 << ", boost factor: " << boostFactor_
1693 << ", algorithm: " << chebyshevAlgorithm_;
1694 if (!userInvDiag_.is_null())
1695 oss << ", diagonal: user-supplied";
1696 oss << "}";
1697 return oss.str();
1698}
1699
1700template<class ScalarType, class MV>
1701void
1703describe (Teuchos::FancyOStream& out,
1704 const Teuchos::EVerbosityLevel verbLevel) const
1705{
1706 using Teuchos::TypeNameTraits;
1707 using std::endl;
1708
1709 const Teuchos::EVerbosityLevel vl =
1710 (verbLevel == Teuchos::VERB_DEFAULT) ? Teuchos::VERB_LOW : verbLevel;
1711 if (vl == Teuchos::VERB_NONE) {
1712 return; // print NOTHING
1713 }
1714
1715 // By convention, describe() starts with a tab.
1716 //
1717 // This does affect all processes on which it's valid to print to
1718 // 'out'. However, it does not actually print spaces to 'out'
1719 // unless operator<< gets called, so it's safe to use on all
1720 // processes.
1721 Teuchos::OSTab tab0 (out);
1722
1723 // We only print on Process 0 of the matrix's communicator. If
1724 // the matrix isn't set, we don't have a communicator, so we have
1725 // to assume that every process can print.
1726 int myRank = -1;
1727 if (A_.is_null () || A_->getComm ().is_null ()) {
1728 myRank = 0;
1729 }
1730 else {
1731 myRank = A_->getComm ()->getRank ();
1732 }
1733 if (myRank == 0) {
1734 // YAML requires quoting the key in this case, to distinguish
1735 // key's colons from the colon that separates key from value.
1736 out << "\"Ifpack2::Details::Chebyshev\":" << endl;
1737 }
1738 Teuchos::OSTab tab1 (out);
1739
1740 if (vl == Teuchos::VERB_LOW) {
1741 if (myRank == 0) {
1742 out << "degree: " << numIters_ << endl
1743 << "lambdaMax: " << lambdaMaxForApply_ << endl
1744 << "alpha: " << eigRatioForApply_ << endl
1745 << "lambdaMin: " << lambdaMinForApply_ << endl
1746 << "boost factor: " << boostFactor_ << endl;
1747 }
1748 return;
1749 }
1750
1751 // vl > Teuchos::VERB_LOW
1752
1753 if (myRank == 0) {
1754 out << "Template parameters:" << endl;
1755 {
1756 Teuchos::OSTab tab2 (out);
1757 out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1758 << "MV: " << TypeNameTraits<MV>::name () << endl;
1759 }
1760
1761 // "Computed parameters" literally means "parameters whose
1762 // values were computed by compute()."
1763 if (myRank == 0) {
1764 out << "Computed parameters:" << endl;
1765 }
1766 }
1767
1768 // Print computed parameters
1769 {
1770 Teuchos::OSTab tab2 (out);
1771 // Users might want to see the values in the computed inverse
1772 // diagonal, so we print them out at the highest verbosity.
1773 if (myRank == 0) {
1774 out << "D_: ";
1775 }
1776 if (D_.is_null ()) {
1777 if (myRank == 0) {
1778 out << "unset" << endl;
1779 }
1780 }
1781 else if (vl <= Teuchos::VERB_HIGH) {
1782 if (myRank == 0) {
1783 out << "set" << endl;
1784 }
1785 }
1786 else { // D_ not null and vl > Teuchos::VERB_HIGH
1787 if (myRank == 0) {
1788 out << endl;
1789 }
1790 // By convention, describe() first indents, then prints.
1791 // We can rely on other describe() implementations to do that.
1792 D_->describe (out, vl);
1793 }
1794 if (myRank == 0) {
1795 // W_ is scratch space; its values are irrelevant.
1796 // All that matters is whether or not they have been set.
1797 out << "W_: " << (W_.is_null () ? "unset" : "set") << endl
1798 << "computedLambdaMax_: " << computedLambdaMax_ << endl
1799 << "computedLambdaMin_: " << computedLambdaMin_ << endl
1800 << "lambdaMaxForApply_: " << lambdaMaxForApply_ << endl
1801 << "lambdaMinForApply_: " << lambdaMinForApply_ << endl
1802 << "eigRatioForApply_: " << eigRatioForApply_ << endl;
1803 }
1804 } // print computed parameters
1805
1806 if (myRank == 0) {
1807 out << "User parameters:" << endl;
1808 }
1809
1810 // Print user parameters
1811 {
1812 Teuchos::OSTab tab2 (out);
1813 out << "userInvDiag_: ";
1814 if (userInvDiag_.is_null ()) {
1815 out << "unset" << endl;
1816 }
1817 else if (vl <= Teuchos::VERB_HIGH) {
1818 out << "set" << endl;
1819 }
1820 else { // userInvDiag_ not null and vl > Teuchos::VERB_HIGH
1821 if (myRank == 0) {
1822 out << endl;
1823 }
1824 userInvDiag_->describe (out, vl);
1825 }
1826 if (myRank == 0) {
1827 out << "userLambdaMax_: " << userLambdaMax_ << endl
1828 << "userLambdaMin_: " << userLambdaMin_ << endl
1829 << "userEigRatio_: " << userEigRatio_ << endl
1830 << "numIters_: " << numIters_ << endl
1831 << "eigMaxIters_: " << eigMaxIters_ << endl
1832 << "eigRelTolerance_: " << eigRelTolerance_ << endl
1833 << "eigNormalizationFreq_: " << eigNormalizationFreq_ << endl
1834 << "zeroStartingSolution_: " << zeroStartingSolution_ << endl
1835 << "assumeMatrixUnchanged_: " << assumeMatrixUnchanged_ << endl
1836 << "chebyshevAlgorithm_: " << chebyshevAlgorithm_ << endl
1837 << "computeMaxResNorm_: " << computeMaxResNorm_ << endl;
1838 }
1839 } // print user parameters
1840}
1841
1842} // namespace Details
1843} // namespace Ifpack2
1844
1845#define IFPACK2_DETAILS_CHEBYSHEV_INSTANT(S,LO,GO,N) \
1846 template class Ifpack2::Details::Chebyshev< S, Tpetra::MultiVector<S, LO, GO, N> >;
1847
1848#endif // IFPACK2_DETAILS_CHEBYSHEV_DEF_HPP
Definition of Chebyshev implementation.
std::vector< ScalarType > optimalWeightsImpl(const int chebyOrder)
Generate optimal weights for using the fourth kind Chebyshev polynomials see: https://arxiv....
Definition Ifpack2_Details_Chebyshev_Weights.hpp:41
Declaration of Chebyshev implementation.
Definition of power methods.
V::scalar_type powerMethodWithInitGuess(const OperatorType &A, const V &D_inv, const int numIters, Teuchos::RCP< V > x, Teuchos::RCP< V > y, const typename Teuchos::ScalarTraits< typename V::scalar_type >::magnitudeType tolerance=1e-7, const int eigNormalizationFreq=1, Teuchos::RCP< Teuchos::FancyOStream > out=Teuchos::null, const bool computeSpectralRadius=true)
Use the power method to estimate the maximum eigenvalue of A*D_inv, given an initial guess vector x.
Definition Ifpack2_PowerMethod.hpp:89
void computeInitialGuessForPowerMethod(V &x, const bool nonnegativeRealParts)
Fill x with random initial guess for power method.
Definition Ifpack2_PowerMethod.hpp:283
Diagonally scaled Chebyshev iteration for Tpetra sparse matrices.
Definition Ifpack2_Chebyshev_decl.hpp:174
MatrixType::scalar_type getLambdaMaxForApply() const
The estimate of the maximum eigenvalue used in the apply().
Definition Ifpack2_Chebyshev_def.hpp:477
Compute scaled damped residual for Chebyshev.
Definition Ifpack2_Details_ChebyshevKernel_decl.hpp:45
Left-scaled Chebyshev iteration.
Definition Ifpack2_Details_Chebyshev_decl.hpp:75
Tpetra::Vector< typename MV::scalar_type, typename MV::local_ordinal_type, typename MV::global_ordinal_type, typename MV::node_type > V
Definition Ifpack2_Details_Chebyshev_decl.hpp:100
Teuchos::ScalarTraits< ScalarType > STS
Traits class for ST.
Definition Ifpack2_Details_Chebyshev_decl.hpp:83
void compute()
(Re)compute the left scaling D_inv, and estimate min and max eigenvalues of D_inv * A.
Definition Ifpack2_Details_Chebyshev_def.hpp:817
void setParameters(Teuchos::ParameterList &plist)
Definition Ifpack2_Details_Chebyshev_def.hpp:349
void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Set the matrix.
Definition Ifpack2_Details_Chebyshev_def.hpp:776
std::string description() const
A single-line description of the Chebyshev solver.
Definition Ifpack2_Details_Chebyshev_def.hpp:1682
Chebyshev(Teuchos::RCP< const row_matrix_type > A)
Definition Ifpack2_Details_Chebyshev_def.hpp:283
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print a description of the Chebyshev solver to out.
Definition Ifpack2_Details_Chebyshev_def.hpp:1703
void print(std::ostream &out)
Print instance data to the given output stream.
Definition Ifpack2_Details_Chebyshev_def.hpp:1066
scalar_type ST
Definition Ifpack2_Details_Chebyshev_decl.hpp:81
STS::magnitudeType MT
Definition Ifpack2_Details_Chebyshev_decl.hpp:85
MT apply(const MV &B, MV &X)
Solve Ax=b for x with Chebyshev iteration with left diagonal scaling.
Definition Ifpack2_Details_Chebyshev_def.hpp:1010
bool hasTransposeApply() const
Whether it's possible to apply the transpose of this operator.
Definition Ifpack2_Details_Chebyshev_def.hpp:1641
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the matrix given to the constructor.
Definition Ifpack2_Details_Chebyshev_def.hpp:1634
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41