Belos Version of the Day
Loading...
Searching...
No Matches
BelosGCRODRSolMgr.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Belos: Block Linear Solvers Package
4//
5// Copyright 2004 NTESS and the Belos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef BELOS_GCRODR_SOLMGR_HPP
11#define BELOS_GCRODR_SOLMGR_HPP
12
16
17#include "BelosConfigDefs.hpp"
21#include "BelosTypes.hpp"
22
23#include "BelosGCRODRIter.hpp"
30#include "Teuchos_BLAS.hpp" // includes Teuchos_ConfigDefs.hpp
31#include "Teuchos_LAPACK.hpp"
32#include "Teuchos_as.hpp"
33#ifdef BELOS_TEUCHOS_TIME_MONITOR
34# include "Teuchos_TimeMonitor.hpp"
35#endif // BELOS_TEUCHOS_TIME_MONITOR
36#if defined(HAVE_TEUCHOSCORE_CXX11)
37# include <type_traits>
38#endif // defined(HAVE_TEUCHOSCORE_CXX11)
39
49
50namespace Belos {
51
53
54
62 public:
63 GCRODRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg) {}
64 };
65
73 public:
74 GCRODRSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg) {}
75 };
76
84 public:
85 GCRODRSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg) {}
86 };
87
95 public:
96 GCRODRSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg) {}
97 };
98
100
124
125 template<class ScalarType, class MV, class OP,
126 const bool lapackSupportsScalarType =
129 public Details::SolverManagerRequiresLapack<ScalarType,MV,OP>
130 {
131 static const bool requiresLapack =
133 typedef Details::SolverManagerRequiresLapack<ScalarType, MV, OP,
134 requiresLapack> base_type;
135
136 public:
138 base_type ()
139 {}
140 GCRODRSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >& problem,
141 const Teuchos::RCP<Teuchos::ParameterList>& pl) :
142 base_type ()
143 {}
144 virtual ~GCRODRSolMgr () {}
145 };
146
150 template<class ScalarType, class MV, class OP>
151 class GCRODRSolMgr<ScalarType, MV, OP, true> :
152 public Details::SolverManagerRequiresLapack<ScalarType, MV, OP, true>
153 {
154
155#if defined(HAVE_TEUCHOSCORE_CXX11)
156# if defined(HAVE_TEUCHOS_COMPLEX)
157 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
158 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
159 std::is_same<ScalarType, std::complex<double> >::value ||
160 std::is_same<ScalarType, Kokkos::complex<double> >::value ||
161 std::is_same<ScalarType, float>::value ||
162 std::is_same<ScalarType, double>::value ||
163 std::is_same<ScalarType, long double>::value,
164 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
165 "types (S,D,C,Z) supported by LAPACK or long double (largely not impl'd).");
166 #else
167 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
168 std::is_same<ScalarType, std::complex<double> >::value ||
169 std::is_same<ScalarType, Kokkos::complex<double> >::value ||
170 std::is_same<ScalarType, float>::value ||
171 std::is_same<ScalarType, double>::value,
172 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
173 "types (S,D,C,Z) supported by LAPACK.");
174 #endif
175# else
176 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
177 static_assert (std::is_same<ScalarType, float>::value ||
178 std::is_same<ScalarType, double>::value ||
179 std::is_same<ScalarType, long double>::value,
180 "Belos::GCRODRSolMgr: ScalarType must be float, double or long double. "
181 "Complex arithmetic support is currently disabled. To "
182 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
183 #else
184 static_assert (std::is_same<ScalarType, float>::value ||
185 std::is_same<ScalarType, double>::value,
186 "Belos::GCRODRSolMgr: ScalarType must be float or double. "
187 "Complex arithmetic support is currently disabled. To "
188 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
189 #endif
190# endif // defined(HAVE_TEUCHOS_COMPLEX)
191#endif // defined(HAVE_TEUCHOSCORE_CXX11)
192
193 private:
196 typedef Teuchos::ScalarTraits<ScalarType> SCT;
197 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
198 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
199 typedef OrthoManagerFactory<ScalarType, MV, OP> ortho_factory_type;
200
201 public:
203
204
210 GCRODRSolMgr();
211
264 GCRODRSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
265 const Teuchos::RCP<Teuchos::ParameterList> &pl);
266
268 virtual ~GCRODRSolMgr() {};
269
271 Teuchos::RCP<SolverManager<ScalarType, MV, OP> > clone () const override {
272 return Teuchos::rcp(new GCRODRSolMgr<ScalarType,MV,OP,true>);
273 }
274
275
277
278
282 return *problem_;
283 }
284
287 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override;
288
291 Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const override {
292 return params_;
293 }
294
300 Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
301 return Teuchos::tuple(timerSolve_);
302 }
303
309 MagnitudeType achievedTol() const override {
310 return achievedTol_;
311 }
312
314 int getNumIters() const override {
315 return numIters_;
316 }
317
320 bool isLOADetected() const override { return false; }
321
323
325
326
328 void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override {
329 problem_ = problem;
330 }
331
333 void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
334
336
338
339
343 void reset( const ResetType type ) override {
344 if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) {
345 bool set = problem_->setProblem();
346 if (!set)
347 throw "Could not set problem.";
348 }
349 else if (type & Belos::RecycleSubspace) {
350 keff = 0;
351 }
352 }
353
354
356
357
384 ReturnType solve() override;
385
387
389
391 std::string description() const override;
392
394
395 private:
396
397 // Called by all constructors; Contains init instructions common to all constructors
398 void init();
399
400 // Initialize solver state storage
401 void initializeStateStorage();
402
403 // Compute updated recycle space given existing recycle space and newly generated Krylov space
404 void buildRecycleSpace2(Teuchos::RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter);
405
406 // Computes harmonic eigenpairs of projected matrix created during the priming solve.
407 // HH is the projected problem from the initial cycle of Gmres, it is (at least) of dimension m+1 x m.
408 // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
409 // The return value is the number of vectors needed to be stored, recycledBlocks or recycledBlocks+1.
410 int getHarmonicVecs1(int m,
411 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
412 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
413
414 // Computes harmonic eigenpairs of projected matrix created during one cycle.
415 // HH is the total block projected problem from the GCRO-DR algorithm, it is (at least) of dimension keff+m+1 x keff+m.
416 // VV is the Krylov vectors from the projected GMRES algorithm, which has (at least) m+1 vectors.
417 // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
418 // The return value is the number of vectors needed to be stored, recycledBlocks or recycledBlocks+1.
419 int getHarmonicVecs2(int keff, int m,
420 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
421 const Teuchos::RCP<const MV>& VV,
422 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
423
424 // Sort list of n floating-point numbers and return permutation vector
425 void sort(std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm);
426
427 // Lapack interface
428 Teuchos::LAPACK<int,ScalarType> lapack;
429
430 // Linear problem.
431 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
432
433 // Output manager.
434 Teuchos::RCP<OutputManager<ScalarType> > printer_;
435 Teuchos::RCP<std::ostream> outputStream_;
436
437 // Status test.
438 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
439 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
440 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
441 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
442 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
443
447 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
448
449 // Current parameter list.
450 Teuchos::RCP<Teuchos::ParameterList> params_;
451
452 // Default solver values.
453 static constexpr double orthoKappa_default_ = 0.0;
454 static constexpr int maxRestarts_default_ = 100;
455 static constexpr int maxIters_default_ = 1000;
456 static constexpr int numBlocks_default_ = 50;
457 static constexpr int blockSize_default_ = 1;
458 static constexpr int recycledBlocks_default_ = 5;
459 static constexpr int verbosity_default_ = Belos::Errors;
460 static constexpr int outputStyle_default_ = Belos::General;
461 static constexpr int outputFreq_default_ = -1;
462 static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
463 static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
464 static constexpr const char * label_default_ = "Belos";
465 static constexpr const char * orthoType_default_ = "ICGS";
466
467 // Current solver values.
468 MagnitudeType convTol_, orthoKappa_, achievedTol_;
469 int maxRestarts_, maxIters_, numIters_;
470 int verbosity_, outputStyle_, outputFreq_;
471 std::string orthoType_;
472 std::string impResScale_, expResScale_;
473
475 // Solver State Storage
477 //
478 // The number of blocks and recycle blocks (m and k, respectively)
479 int numBlocks_, recycledBlocks_;
480 // Current size of recycled subspace
481 int keff;
482 //
483 // Residual vector
484 Teuchos::RCP<MV> r_;
485 //
486 // Search space
487 Teuchos::RCP<MV> V_;
488 //
489 // Recycled subspace and its image
490 Teuchos::RCP<MV> U_, C_;
491 //
492 // Updated recycle space and its image
493 Teuchos::RCP<MV> U1_, C1_;
494 //
495 // Storage used in constructing harmonic Ritz values/vectors
496 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H2_;
497 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
498 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B_;
499 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PP_;
500 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > HP_;
501 std::vector<ScalarType> tau_;
502 std::vector<ScalarType> work_;
503 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
504 std::vector<int> ipiv_;
506
507 // Timers.
508 std::string label_;
509 Teuchos::RCP<Teuchos::Time> timerSolve_;
510
511 // Internal state variables.
512 bool isSet_;
513
514 // Have we generated or regenerated a recycle space yet this solve?
515 bool builtRecycleSpace_;
516 };
517
518
519// Empty Constructor
520template<class ScalarType, class MV, class OP>
522 achievedTol_(0.0),
523 numIters_(0)
524{
525 init ();
526}
527
528
529// Basic Constructor
530template<class ScalarType, class MV, class OP>
532GCRODRSolMgr(const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >& problem,
533 const Teuchos::RCP<Teuchos::ParameterList>& pl):
534 achievedTol_(0.0),
535 numIters_(0)
536{
537 // Initialize local pointers to null, and initialize local variables
538 // to default values.
539 init ();
540
541 TEUCHOS_TEST_FOR_EXCEPTION(
542 problem == Teuchos::null, std::invalid_argument,
543 "Belos::GCRODRSolMgr constructor: The solver manager's "
544 "constructor needs the linear problem argument 'problem' "
545 "to be non-null.");
546 problem_ = problem;
547
548 // Set the parameters using the list that was passed in. If null,
549 // we defer initialization until a non-null list is set (by the
550 // client calling setParameters(), or by calling solve() -- in
551 // either case, a null parameter list indicates that default
552 // parameters should be used).
553 if (! pl.is_null ()) {
554 setParameters (pl);
555 }
556}
557
558// Common instructions executed in all constructors
559template<class ScalarType, class MV, class OP>
561 outputStream_ = Teuchos::rcpFromRef(std::cout);
563 orthoKappa_ = orthoKappa_default_;
564 maxRestarts_ = maxRestarts_default_;
565 maxIters_ = maxIters_default_;
566 numBlocks_ = numBlocks_default_;
567 recycledBlocks_ = recycledBlocks_default_;
568 verbosity_ = verbosity_default_;
569 outputStyle_ = outputStyle_default_;
570 outputFreq_ = outputFreq_default_;
571 orthoType_ = orthoType_default_;
572 impResScale_ = impResScale_default_;
573 expResScale_ = expResScale_default_;
574 label_ = label_default_;
575 isSet_ = false;
576 builtRecycleSpace_ = false;
577 keff = 0;
578 r_ = Teuchos::null;
579 V_ = Teuchos::null;
580 U_ = Teuchos::null;
581 C_ = Teuchos::null;
582 U1_ = Teuchos::null;
583 C1_ = Teuchos::null;
584 PP_ = Teuchos::null;
585 HP_ = Teuchos::null;
586 H2_ = Teuchos::null;
587 R_ = Teuchos::null;
588 H_ = Teuchos::null;
589 B_ = Teuchos::null;
590}
591
592template<class ScalarType, class MV, class OP>
593void
595setParameters (const Teuchos::RCP<Teuchos::ParameterList> &params)
596{
597 using Teuchos::isParameterType;
598 using Teuchos::getParameter;
599 using Teuchos::null;
600 using Teuchos::ParameterList;
601 using Teuchos::parameterList;
602 using Teuchos::RCP;
603 using Teuchos::rcp;
604 using Teuchos::rcp_dynamic_cast;
605 using Teuchos::rcpFromRef;
606 using Teuchos::Exceptions::InvalidParameter;
607 using Teuchos::Exceptions::InvalidParameterName;
608 using Teuchos::Exceptions::InvalidParameterType;
609
610 // The default parameter list contains all parameters that
611 // GCRODRSolMgr understands, and none that it doesn't understand.
612 RCP<const ParameterList> defaultParams = getValidParameters();
613
614 // Create the internal parameter list if one doesn't already exist.
615 //
616 // (mfh 28 Feb 2011, 10 Mar 2011) At the time this code was written,
617 // ParameterList did not have validators or validateParameters().
618 // This is why the code below carefully validates the parameters one
619 // by one and fills in defaults. This code could be made a lot
620 // shorter by using validators. To do so, we would have to define
621 // appropriate validators for all the parameters. (This would more
622 // or less just move all that validation code out of this routine
623 // into to getValidParameters().)
624 //
625 // For an analogous reason, GCRODRSolMgr defines default parameter
626 // values as class data, as well as in the default ParameterList.
627 // This redundancy could be removed by defining the default
628 // parameter values only in the default ParameterList (which
629 // documents each parameter as well -- handy!).
630 if (params_.is_null()) {
631 params_ = parameterList (*defaultParams);
632 } else {
633 // A common case for setParameters() is for it to be called at the
634 // beginning of the solve() routine. This follows the Belos
635 // pattern of delaying initialization until the last possible
636 // moment (when the user asks Belos to perform the solve). In
637 // this common case, we save ourselves a deep copy of the input
638 // parameter list.
639 if (params_ != params) {
640 // Make a deep copy of the input parameter list. This allows
641 // the caller to modify or change params later, without
642 // affecting the behavior of this solver. This solver will then
643 // only change its internal parameters if setParameters() is
644 // called again.
645 params_ = parameterList (*params);
646 }
647
648 // Fill in any missing parameters and their default values. Also,
649 // throw an exception if the parameter list has any misspelled or
650 // "extra" parameters. If you don't like this behavior, you'll
651 // want to replace the line of code below with your desired
652 // validation scheme. Note that Teuchos currently only implements
653 // two options:
654 //
655 // 1. validateParameters() requires that params_ has all the
656 // parameters that the default list has, and none that it
657 // doesn't have.
658 //
659 // 2. validateParametersAndSetDefaults() fills in missing
660 // parameters in params_ using the default list, but requires
661 // that any parameter provided in params_ is also in the
662 // default list.
663 //
664 // Here is an easy way to ignore any "extra" or misspelled
665 // parameters: Make a deep copy of the default list, fill in any
666 // "missing" parameters from the _input_ list, and then validate
667 // the input list using the deep copy of the default list. We
668 // show this in code:
669 //
670 // RCP<ParameterList> defaultCopy = parameterList (*getValidParameters ());
671 // defaultCopy->validateParametersAndSetDefaults (params);
672 // params->validateParametersAndSetDefaults (defaultCopy);
673 //
674 // This method is not entirely robust, because the input list may
675 // have incorrect validators set for existing parameters in the
676 // default list. This would then cause "validation" of the
677 // default list to throw an exception. As a result, we've chosen
678 // for now to be intolerant of misspellings and "extra" parameters
679 // in the input list.
680 params_->validateParametersAndSetDefaults (*defaultParams);
681 }
682
683 // Check for maximum number of restarts.
684 if (params->isParameter ("Maximum Restarts")) {
685 maxRestarts_ = params->get("Maximum Restarts", maxRestarts_default_);
686
687 // Update parameter in our list.
688 params_->set ("Maximum Restarts", maxRestarts_);
689 }
690
691 // Check for maximum number of iterations
692 if (params->isParameter ("Maximum Iterations")) {
693 maxIters_ = params->get ("Maximum Iterations", maxIters_default_);
694
695 // Update parameter in our list and in status test.
696 params_->set ("Maximum Iterations", maxIters_);
697 if (! maxIterTest_.is_null())
698 maxIterTest_->setMaxIters (maxIters_);
699 }
700
701 // Check for the maximum number of blocks.
702 if (params->isParameter ("Num Blocks")) {
703 numBlocks_ = params->get ("Num Blocks", numBlocks_default_);
704 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
705 "Belos::GCRODRSolMgr: The \"Num Blocks\" parameter must "
706 "be strictly positive, but you specified a value of "
707 << numBlocks_ << ".");
708 // Update parameter in our list.
709 params_->set ("Num Blocks", numBlocks_);
710 }
711
712 // Check for the maximum number of blocks.
713 if (params->isParameter ("Num Recycled Blocks")) {
714 recycledBlocks_ = params->get ("Num Recycled Blocks",
715 recycledBlocks_default_);
716 TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks_ <= 0, std::invalid_argument,
717 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
718 "parameter must be strictly positive, but you specified "
719 "a value of " << recycledBlocks_ << ".");
720 TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks_ >= numBlocks_, std::invalid_argument,
721 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
722 "parameter must be less than the \"Num Blocks\" "
723 "parameter, but you specified \"Num Recycled Blocks\" "
724 "= " << recycledBlocks_ << " and \"Num Blocks\" = "
725 << numBlocks_ << ".");
726 // Update parameter in our list.
727 params_->set("Num Recycled Blocks", recycledBlocks_);
728 }
729
730 // Check to see if the timer label changed. If it did, update it in
731 // the parameter list, and create a new timer with that label (if
732 // Belos was compiled with timers enabled).
733 if (params->isParameter ("Timer Label")) {
734 std::string tempLabel = params->get ("Timer Label", label_default_);
735
736 // Update parameter in our list and solver timer
737 if (tempLabel != label_) {
738 label_ = tempLabel;
739 params_->set ("Timer Label", label_);
740 std::string solveLabel = label_ + ": GCRODRSolMgr total solve time";
741#ifdef BELOS_TEUCHOS_TIME_MONITOR
742 timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
743#endif
744 if (ortho_ != Teuchos::null) {
745 ortho_->setLabel( label_ );
746 }
747 }
748 }
749
750 // Check for a change in verbosity level
751 if (params->isParameter ("Verbosity")) {
752 if (isParameterType<int> (*params, "Verbosity")) {
753 verbosity_ = params->get ("Verbosity", verbosity_default_);
754 } else {
755 verbosity_ = (int) getParameter<Belos::MsgType> (*params, "Verbosity");
756 }
757 // Update parameter in our list.
758 params_->set ("Verbosity", verbosity_);
759 // If the output manager (printer_) is null, then we will
760 // instantiate it later with the correct verbosity.
761 if (! printer_.is_null())
762 printer_->setVerbosity (verbosity_);
763 }
764
765 // Check for a change in output style
766 if (params->isParameter ("Output Style")) {
767 if (isParameterType<int> (*params, "Output Style")) {
768 outputStyle_ = params->get ("Output Style", outputStyle_default_);
769 } else {
770 outputStyle_ = (int) getParameter<OutputType> (*params, "Output Style");
771 }
772
773 // Update parameter in our list.
774 params_->set ("Output Style", outputStyle_);
775 // We will (re)instantiate the output status test afresh below.
776 outputTest_ = null;
777 }
778
779 // Get the output stream for the output manager.
780 //
781 // While storing the output stream in the parameter list (either as
782 // an RCP or as a nonconst reference) is convenient and safe for
783 // programming, it makes it impossible to serialize the parameter
784 // list, read it back in from the serialized representation, and get
785 // the same output stream as before. This is because output streams
786 // may be arbitrary constructed objects.
787 //
788 // In case the user tried reading in the parameter list from a
789 // serialized representation and the output stream can't be read
790 // back in, we set the output stream to point to std::cout. This
791 // ensures reasonable behavior.
792 if (params->isParameter ("Output Stream")) {
793 try {
794 outputStream_ = getParameter<RCP<std::ostream> > (*params, "Output Stream");
795 } catch (InvalidParameter&) {
796 outputStream_ = rcpFromRef (std::cout);
797 }
798 // We assume that a null output stream indicates that the user
799 // doesn't want to print anything, so we replace it with a "black
800 // hole" stream that prints nothing sent to it. (We can't use a
801 // null output stream, since the output manager always sends
802 // things it wants to print to the output stream.)
803 if (outputStream_.is_null()) {
804 outputStream_ = rcp (new Teuchos::oblackholestream);
805 }
806 // Update parameter in our list.
807 params_->set ("Output Stream", outputStream_);
808 // If the output manager (printer_) is null, then we will
809 // instantiate it later with the correct output stream.
810 if (! printer_.is_null()) {
811 printer_->setOStream (outputStream_);
812 }
813 }
814
815 // frequency level
816 if (verbosity_ & Belos::StatusTestDetails) {
817 if (params->isParameter ("Output Frequency")) {
818 outputFreq_ = params->get ("Output Frequency", outputFreq_default_);
819 }
820
821 // Update parameter in out list and output status test.
822 params_->set("Output Frequency", outputFreq_);
823 if (! outputTest_.is_null())
824 outputTest_->setOutputFrequency (outputFreq_);
825 }
826
827 // Create output manager if we need to, using the verbosity level
828 // and output stream that we fetched above. We do this here because
829 // instantiating an OrthoManager using OrthoManagerFactory requires
830 // a valid OutputManager.
831 if (printer_.is_null()) {
832 printer_ = rcp (new OutputManager<ScalarType> (verbosity_, outputStream_));
833 }
834
835 // Get the orthogonalization manager name ("Orthogonalization").
836 //
837 // Getting default values for the orthogonalization manager
838 // parameters ("Orthogonalization Parameters") requires knowing the
839 // orthogonalization manager name. Save it for later, and also
840 // record whether it's different than before.
842 bool changedOrthoType = false;
843 if (params->isParameter ("Orthogonalization")) {
844 const std::string& tempOrthoType =
845 params->get ("Orthogonalization", orthoType_default_);
846 // Ensure that the specified orthogonalization type is valid.
847 if (! factory.isValidName (tempOrthoType)) {
848 std::ostringstream os;
849 os << "Belos::GCRODRSolMgr: Invalid orthogonalization name \""
850 << tempOrthoType << "\". The following are valid options "
851 << "for the \"Orthogonalization\" name parameter: ";
852 factory.printValidNames (os);
853 throw std::invalid_argument (os.str());
854 }
855 if (tempOrthoType != orthoType_) {
856 changedOrthoType = true;
857 orthoType_ = tempOrthoType;
858 // Update parameter in our list.
859 params_->set ("Orthogonalization", orthoType_);
860 }
861 }
862
863 // Get any parameters for the orthogonalization ("Orthogonalization
864 // Parameters"). If not supplied, the orthogonalization manager
865 // factory will supply default values.
866 //
867 // NOTE (mfh 12 Jan 2011) For the sake of backwards compatibility,
868 // if params has an "Orthogonalization Constant" parameter and the
869 // DGKS orthogonalization manager is to be used, the value of this
870 // parameter will override DGKS's "depTol" parameter.
871 //
872 // Users must supply the orthogonalization manager parameters as a
873 // sublist (supplying it as an RCP<ParameterList> would make the
874 // resulting parameter list not serializable).
875 RCP<ParameterList> orthoParams;
876 { // The nonmember function sublist() returns an RCP<ParameterList>,
877 // which is what we want here.
878 using Teuchos::sublist;
879 // Abbreviation to avoid typos.
880 const std::string paramName ("Orthogonalization Parameters");
881
882 try {
883 orthoParams = sublist (params_, paramName, true);
884 } catch (InvalidParameter&) {
885 // We didn't get the parameter list from params, so get a
886 // default parameter list from the OrthoManagerFactory. Modify
887 // params_ so that it has the default parameter list, and set
888 // orthoParams to ensure it's a sublist of params_ (and not just
889 // a copy of one).
890 params_->set (paramName, factory.getDefaultParameters (orthoType_));
891 orthoParams = sublist (params_, paramName, true);
892 }
893 }
894 TEUCHOS_TEST_FOR_EXCEPTION(orthoParams.is_null(), std::logic_error,
895 "Failed to get orthogonalization parameters. "
896 "Please report this bug to the Belos developers.");
897
898 // Instantiate a new MatOrthoManager subclass instance if necessary.
899 // If not necessary, then tell the existing instance about the new
900 // parameters.
901 if (ortho_.is_null() || changedOrthoType) {
902 // We definitely need to make a new MatOrthoManager, since either
903 // we haven't made one yet, or we've changed orthogonalization
904 // methods. Creating the orthogonalization manager requires that
905 // the OutputManager (printer_) already be initialized.
906 ortho_ = factory.makeMatOrthoManager (orthoType_, null, printer_,
907 label_, orthoParams);
908 } else {
909 // If the MatOrthoManager implements the ParameterListAcceptor
910 // mix-in interface, we can propagate changes to its parameters
911 // without reinstantiating the MatOrthoManager.
912 //
913 // We recommend that all MatOrthoManager subclasses implement
914 // Teuchos::ParameterListAcceptor, but do not (yet) require this.
915 typedef Teuchos::ParameterListAcceptor PLA;
916 RCP<PLA> pla = rcp_dynamic_cast<PLA> (ortho_);
917 if (pla.is_null()) {
918 // Oops, it's not a ParameterListAcceptor. We have to
919 // reinstantiate the MatOrthoManager in order to pass in the
920 // possibly new parameters.
921 ortho_ = factory.makeMatOrthoManager (orthoType_, null, printer_,
922 label_, orthoParams);
923 } else {
924 pla->setParameterList (orthoParams);
925 }
926 }
927
928 // The DGKS orthogonalization accepts a "Orthogonalization Constant"
929 // parameter (also called kappa in the code, but not in the
930 // parameter list). If its value is provided in the given parameter
931 // list, and its value is positive, use it. Ignore negative values.
932 //
933 // NOTE (mfh 12 Jan 2011) This overrides the "depTol" parameter that
934 // may have been specified in "Orthogonalization Parameters". We
935 // retain this behavior for backwards compatibility.
936 if (params->isParameter ("Orthogonalization Constant")) {
937 MagnitudeType orthoKappa = orthoKappa_default_;
938 if (params->isType<MagnitudeType> ("Orthogonalization Constant")) {
939 orthoKappa = params->get ("Orthogonalization Constant", orthoKappa);
940 }
941 else {
942 orthoKappa = params->get ("Orthogonalization Constant", orthoKappa_default_);
943 }
944
945 if (orthoKappa > 0) {
946 orthoKappa_ = orthoKappa;
947 // Update parameter in our list.
948 params_->set("Orthogonalization Constant", orthoKappa_);
949 // Only DGKS currently accepts this parameter.
950 if (orthoType_ == "DGKS" && ! ortho_.is_null()) {
951 typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_man_type;
952 // This cast should always succeed; it's a bug
953 // otherwise. (If the cast fails, then orthoType_
954 // doesn't correspond to the OrthoManager subclass
955 // instance that we think we have, so we initialized the
956 // wrong subclass somehow.)
957 rcp_dynamic_cast<ortho_man_type>(ortho_)->setDepTol (orthoKappa_);
958 }
959 }
960 }
961
962 // Convergence
963 typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
964 typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
965
966 // Check for convergence tolerance
967 if (params->isParameter("Convergence Tolerance")) {
968 if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
969 convTol_ = params->get ("Convergence Tolerance",
970 static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
971 }
972 else {
973 convTol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
974 }
975
976 // Update parameter in our list and residual tests.
977 params_->set ("Convergence Tolerance", convTol_);
978 if (! impConvTest_.is_null())
979 impConvTest_->setTolerance (convTol_);
980 if (! expConvTest_.is_null())
981 expConvTest_->setTolerance (convTol_);
982 }
983
984 // Check for a change in scaling, if so we need to build new residual tests.
985 if (params->isParameter ("Implicit Residual Scaling")) {
986 std::string tempImpResScale =
987 getParameter<std::string> (*params, "Implicit Residual Scaling");
988
989 // Only update the scaling if it's different.
990 if (impResScale_ != tempImpResScale) {
991 ScaleType impResScaleType = convertStringToScaleType (tempImpResScale);
992 impResScale_ = tempImpResScale;
993
994 // Update parameter in our list and residual tests
995 params_->set("Implicit Residual Scaling", impResScale_);
996 // NOTE (mfh 28 Feb 2011) StatusTestImpResNorm only lets you
997 // call defineScaleForm() once. The code below attempts to call
998 // defineScaleForm(); if the scale form has already been
999 // defined, it constructs a new StatusTestImpResNorm instance.
1000 // StatusTestImpResNorm really should not expose the
1001 // defineScaleForm() method, since it's serving an
1002 // initialization purpose; all initialization should happen in
1003 // the constructor whenever possible. In that case, the code
1004 // below could be simplified into a single (re)instantiation.
1005 if (! impConvTest_.is_null()) {
1006 try {
1007 impConvTest_->defineScaleForm (impResScaleType, Belos::TwoNorm);
1008 }
1009 catch (StatusTestError&) {
1010 // Delete the convergence test so it gets constructed again.
1011 impConvTest_ = null;
1012 convTest_ = null;
1013 }
1014 }
1015 }
1016 }
1017
1018 if (params->isParameter("Explicit Residual Scaling")) {
1019 std::string tempExpResScale =
1020 getParameter<std::string> (*params, "Explicit Residual Scaling");
1021
1022 // Only update the scaling if it's different.
1023 if (expResScale_ != tempExpResScale) {
1024 ScaleType expResScaleType = convertStringToScaleType (tempExpResScale);
1025 expResScale_ = tempExpResScale;
1026
1027 // Update parameter in our list and residual tests
1028 params_->set("Explicit Residual Scaling", expResScale_);
1029 // NOTE (mfh 28 Feb 2011) See note above on the (re)construction
1030 // of StatusTestImpResNorm.
1031 if (! expConvTest_.is_null()) {
1032 try {
1033 expConvTest_->defineScaleForm (expResScaleType, Belos::TwoNorm);
1034 }
1035 catch (StatusTestError&) {
1036 // Delete the convergence test so it gets constructed again.
1037 expConvTest_ = null;
1038 convTest_ = null;
1039 }
1040 }
1041 }
1042 }
1043 //
1044 // Create iteration stopping criteria ("status tests") if we need
1045 // to, by combining three different stopping criteria.
1046 //
1047 // First, construct maximum-number-of-iterations stopping criterion.
1048 if (maxIterTest_.is_null())
1049 maxIterTest_ = rcp (new StatusTestMaxIters<ScalarType,MV,OP> (maxIters_));
1050
1051 // Implicit residual test, using the native residual to determine if
1052 // convergence was achieved.
1053 if (impConvTest_.is_null()) {
1054 impConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
1055 impConvTest_->defineScaleForm (convertStringToScaleType (impResScale_),
1057 }
1058
1059 // Explicit residual test once the native residual is below the tolerance
1060 if (expConvTest_.is_null()) {
1061 expConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
1062 expConvTest_->defineResForm (StatusTestResNorm_t::Explicit, Belos::TwoNorm);
1063 expConvTest_->defineScaleForm (convertStringToScaleType (expResScale_),
1065 }
1066 // Convergence test first tests the implicit residual, then the
1067 // explicit residual if the implicit residual test passes.
1068 if (convTest_.is_null()) {
1069 convTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::SEQ,
1070 impConvTest_,
1071 expConvTest_));
1072 }
1073 // Construct the complete iteration stopping criterion:
1074 //
1075 // "Stop iterating if the maximum number of iterations has been
1076 // reached, or if the convergence test passes."
1077 sTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::OR,
1078 maxIterTest_,
1079 convTest_));
1080 // Create the status test output class.
1081 // This class manages and formats the output from the status test.
1082 StatusTestOutputFactory<ScalarType,MV,OP> stoFactory (outputStyle_);
1083 outputTest_ = stoFactory.create (printer_, sTest_, outputFreq_,
1085
1086 // Set the solver string for the output test
1087 std::string solverDesc = " GCRODR ";
1088 outputTest_->setSolverDesc( solverDesc );
1089
1090 // Create the timer if we need to.
1091 if (timerSolve_.is_null()) {
1092 std::string solveLabel = label_ + ": GCRODRSolMgr total solve time";
1093#ifdef BELOS_TEUCHOS_TIME_MONITOR
1094 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
1095#endif
1096 }
1097
1098 // Inform the solver manager that the current parameters were set.
1099 isSet_ = true;
1100}
1101
1102
1103template<class ScalarType, class MV, class OP>
1104Teuchos::RCP<const Teuchos::ParameterList>
1106{
1107 using Teuchos::ParameterList;
1108 using Teuchos::parameterList;
1109 using Teuchos::RCP;
1110
1111 static RCP<const ParameterList> validPL;
1112 if (is_null(validPL)) {
1113 RCP<ParameterList> pl = parameterList ();
1114
1115 // Set all the valid parameters and their default values.
1116 pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
1117 "The relative residual tolerance that needs to be achieved by the\n"
1118 "iterative solver in order for the linear system to be declared converged.");
1119 pl->set("Maximum Restarts", static_cast<int>(maxRestarts_default_),
1120 "The maximum number of cycles allowed for each\n"
1121 "set of RHS solved.");
1122 pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
1123 "The maximum number of iterations allowed for each\n"
1124 "set of RHS solved.");
1125 // mfh 25 Oct 2010: "Block Size" must be 1 because GCRODR is
1126 // currently not a block method: i.e., it does not work on
1127 // multiple right-hand sides at once.
1128 pl->set("Block Size", static_cast<int>(blockSize_default_),
1129 "Block Size Parameter -- currently must be 1 for GCRODR");
1130 pl->set("Num Blocks", static_cast<int>(numBlocks_default_),
1131 "The maximum number of vectors allowed in the Krylov subspace\n"
1132 "for each set of RHS solved.");
1133 pl->set("Num Recycled Blocks", static_cast<int>(recycledBlocks_default_),
1134 "The maximum number of vectors in the recycled subspace." );
1135 pl->set("Verbosity", static_cast<int>(verbosity_default_),
1136 "What type(s) of solver information should be outputted\n"
1137 "to the output stream.");
1138 pl->set("Output Style", static_cast<int>(outputStyle_default_),
1139 "What style is used for the solver information outputted\n"
1140 "to the output stream.");
1141 pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
1142 "How often convergence information should be outputted\n"
1143 "to the output stream.");
1144 pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
1145 "A reference-counted pointer to the output stream where all\n"
1146 "solver output is sent.");
1147 pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
1148 "The type of scaling used in the implicit residual convergence test.");
1149 pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
1150 "The type of scaling used in the explicit residual convergence test.");
1151 pl->set("Timer Label", static_cast<const char *>(label_default_),
1152 "The string to use as a prefix for the timer labels.");
1153 {
1155 pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
1156 "The type of orthogonalization to use. Valid options: " +
1157 factory.validNamesString());
1158 RCP<const ParameterList> orthoParams =
1159 factory.getDefaultParameters (orthoType_default_);
1160 pl->set ("Orthogonalization Parameters", *orthoParams,
1161 "Parameters specific to the type of orthogonalization used.");
1162 }
1163 pl->set("Orthogonalization Constant",static_cast<MagnitudeType>(orthoKappa_default_),
1164 "When using DGKS orthogonalization: the \"depTol\" constant, used "
1165 "to determine whether another step of classical Gram-Schmidt is "
1166 "necessary. Otherwise ignored.");
1167 validPL = pl;
1168 }
1169 return validPL;
1170}
1171
1172// initializeStateStorage
1173template<class ScalarType, class MV, class OP>
1175
1176 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1177
1178 // Check if there is any multivector to clone from.
1179 Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1180 if (rhsMV == Teuchos::null) {
1181 // Nothing to do
1182 return;
1183 }
1184 else {
1185
1186 // Initialize the state storage
1187 TEUCHOS_TEST_FOR_EXCEPTION(static_cast<ptrdiff_t>(numBlocks_) > MVT::GetGlobalLength(*rhsMV),std::invalid_argument,
1188 "Belos::GCRODRSolMgr::initializeStateStorage(): Cannot generate a Krylov basis with dimension larger the operator!");
1189
1190 // If the subspace has not been initialized before, generate it using the RHS from lp_.
1191 if (U_ == Teuchos::null) {
1192 U_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1193 }
1194 else {
1195 // Generate U_ by cloning itself ONLY if more space is needed.
1196 if (MVT::GetNumberVecs(*U_) < recycledBlocks_+1) {
1197 Teuchos::RCP<const MV> tmp = U_;
1198 U_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1199 }
1200 }
1201
1202 // If the subspace has not been initialized before, generate it using the RHS from lp_.
1203 if (C_ == Teuchos::null) {
1204 C_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1205 }
1206 else {
1207 // Generate C_ by cloning itself ONLY if more space is needed.
1208 if (MVT::GetNumberVecs(*C_) < recycledBlocks_+1) {
1209 Teuchos::RCP<const MV> tmp = C_;
1210 C_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1211 }
1212 }
1213
1214 // If the subspace has not been initialized before, generate it using the RHS from lp_.
1215 if (V_ == Teuchos::null) {
1216 V_ = MVT::Clone( *rhsMV, numBlocks_+1 );
1217 }
1218 else {
1219 // Generate V_ by cloning itself ONLY if more space is needed.
1220 if (MVT::GetNumberVecs(*V_) < numBlocks_+1) {
1221 Teuchos::RCP<const MV> tmp = V_;
1222 V_ = MVT::Clone( *tmp, numBlocks_+1 );
1223 }
1224 }
1225
1226 // If the subspace has not been initialized before, generate it using the RHS from lp_.
1227 if (U1_ == Teuchos::null) {
1228 U1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1229 }
1230 else {
1231 // Generate U1_ by cloning itself ONLY if more space is needed.
1232 if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1233 Teuchos::RCP<const MV> tmp = U1_;
1234 U1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1235 }
1236 }
1237
1238 // If the subspace has not been initialized before, generate it using the RHS from lp_.
1239 if (C1_ == Teuchos::null) {
1240 C1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1241 }
1242 else {
1243 // Generate C1_ by cloning itself ONLY if more space is needed.
1244 if (MVT::GetNumberVecs(*C1_) < recycledBlocks_+1) {
1245 Teuchos::RCP<const MV> tmp = C1_;
1246 C1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1247 }
1248 }
1249
1250 // Generate r_ only if it doesn't exist
1251 if (r_ == Teuchos::null)
1252 r_ = MVT::Clone( *rhsMV, 1 );
1253
1254 // Size of tau_ will change during computation, so just be sure it starts with appropriate size
1255 tau_.resize(recycledBlocks_+1);
1256
1257 // Size of work_ will change during computation, so just be sure it starts with appropriate size
1258 work_.resize(recycledBlocks_+1);
1259
1260 // Size of ipiv_ will change during computation, so just be sure it starts with appropriate size
1261 ipiv_.resize(recycledBlocks_+1);
1262
1263 // Generate H2_ only if it doesn't exist, otherwise resize it.
1264 if (H2_ == Teuchos::null)
1265 H2_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 ) );
1266 else {
1267 if ( (H2_->numRows() != numBlocks_+recycledBlocks_+2) || (H2_->numCols() != numBlocks_+recycledBlocks_+1) )
1268 H2_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1269 }
1270 H2_->putScalar(zero);
1271
1272 // Generate R_ only if it doesn't exist, otherwise resize it.
1273 if (R_ == Teuchos::null)
1274 R_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycledBlocks_+1, recycledBlocks_+1 ) );
1275 else {
1276 if ( (R_->numRows() != recycledBlocks_+1) || (R_->numCols() != recycledBlocks_+1) )
1277 R_->reshape( recycledBlocks_+1, recycledBlocks_+1 );
1278 }
1279 R_->putScalar(zero);
1280
1281 // Generate PP_ only if it doesn't exist, otherwise resize it.
1282 if (PP_ == Teuchos::null)
1283 PP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, recycledBlocks_+1 ) );
1284 else {
1285 if ( (PP_->numRows() != numBlocks_+recycledBlocks_+2) || (PP_->numCols() != recycledBlocks_+1) )
1286 PP_->reshape( numBlocks_+recycledBlocks_+2, recycledBlocks_+1 );
1287 }
1288
1289 // Generate HP_ only if it doesn't exist, otherwise resize it.
1290 if (HP_ == Teuchos::null)
1291 HP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 ) );
1292 else {
1293 if ( (HP_->numRows() != numBlocks_+recycledBlocks_+2) || (HP_->numCols() != numBlocks_+recycledBlocks_+1) )
1294 HP_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1295 }
1296
1297 } // end else
1298}
1299
1300
1301// solve()
1302template<class ScalarType, class MV, class OP>
1304 using Teuchos::RCP;
1305 using Teuchos::rcp;
1306
1307 // Set the current parameters if they were not set before.
1308 // NOTE: This may occur if the user generated the solver manager with the default constructor and
1309 // then didn't set any parameters using setParameters().
1310 if (!isSet_) { setParameters( params_ ); }
1311
1312 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1313 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1314 std::vector<int> index(numBlocks_+1);
1315
1316 TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,GCRODRSolMgrLinearProblemFailure, "Belos::GCRODRSolMgr::solve(): Linear problem is not a valid object.");
1317
1318 TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(),GCRODRSolMgrLinearProblemFailure,"Belos::GCRODRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
1319
1320 // Create indices for the linear systems to be solved.
1321 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1322 std::vector<int> currIdx(1);
1323 currIdx[0] = 0;
1324
1325 // Inform the linear problem of the current linear system to solve.
1326 problem_->setLSIndex( currIdx );
1327
1328 // Check the number of blocks and change them is necessary.
1329 ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1330 if (static_cast<ptrdiff_t>(numBlocks_) > dim) {
1331 numBlocks_ = Teuchos::as<int>(dim);
1332 printer_->stream(Warnings) <<
1333 "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
1334 " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << numBlocks_ << std::endl;
1335 params_->set("Num Blocks", numBlocks_);
1336 }
1337
1338 // Assume convergence is achieved, then let any failed convergence set this to false.
1339 bool isConverged = true;
1340
1341 // Initialize storage for all state variables
1342 initializeStateStorage();
1343
1345 // Parameter list
1346 Teuchos::ParameterList plist;
1347
1348 plist.set("Num Blocks",numBlocks_);
1349 plist.set("Recycled Blocks",recycledBlocks_);
1350
1352 // GCRODR solver
1353
1354 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter;
1355 gcrodr_iter = rcp( new GCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
1356 // Number of iterations required to generate initial recycle space (if needed)
1357 int prime_iterations = 0;
1358
1359 // Enter solve() iterations
1360 {
1361#ifdef BELOS_TEUCHOS_TIME_MONITOR
1362 Teuchos::TimeMonitor slvtimer(*timerSolve_);
1363#endif
1364
1365 while ( numRHS2Solve > 0 ) {
1366
1367 // Set flag indicating recycle space has not been generated this solve
1368 builtRecycleSpace_ = false;
1369
1370 // Reset the status test.
1371 outputTest_->reset();
1372
1374 // Initialize recycled subspace for GCRODR
1375
1376 // If there is a subspace to recycle, recycle it, otherwise generate the initial recycled subspace.
1377 if (keff > 0) {
1378 TEUCHOS_TEST_FOR_EXCEPTION(keff < recycledBlocks_,GCRODRSolMgrRecyclingFailure,
1379 "Belos::GCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
1380
1381 printer_->stream(Debug) << " Now solving RHS index " << currIdx[0] << " using recycled subspace of dimension " << keff << std::endl << std::endl;
1382 // Compute image of U_ under the new operator
1383 index.resize(keff);
1384 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1385 RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1386 RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1387 problem_->apply( *Utmp, *Ctmp );
1388
1389 RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1390
1391 // Orthogonalize this block
1392 // Get a matrix to hold the orthonormalization coefficients.
1393 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff, keff );
1394 int rank = ortho_->normalize(*Ctmp, rcp(&Rtmp,false));
1395 // Throw an error if we could not orthogonalize this block
1396 TEUCHOS_TEST_FOR_EXCEPTION(rank != keff,GCRODRSolMgrOrthoFailure,"Belos::GCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1397
1398 // U_ = U_*R^{-1}
1399 // First, compute LU factorization of R
1400 int info = 0;
1401 ipiv_.resize(Rtmp.numRows());
1402 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1403 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1404
1405 // Now, form inv(R)
1406 int lwork = Rtmp.numRows();
1407 work_.resize(lwork);
1408 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1409 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1410
1411 // U_ = U1_; (via a swap)
1412 MVT::MvTimesMatAddMv( one, *Utmp, Rtmp, zero, *U1tmp );
1413 std::swap(U_, U1_);
1414
1415 // Must reinitialize after swap
1416 index.resize(keff);
1417 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1418 Ctmp = MVT::CloneViewNonConst( *C_, index );
1419 Utmp = MVT::CloneView( *U_, index );
1420
1421 // Compute C_'*r_
1422 Teuchos::SerialDenseMatrix<int,ScalarType> Ctr(keff,1);
1423 problem_->computeCurrPrecResVec( &*r_ );
1424 MVT::MvTransMv( one, *Ctmp, *r_, Ctr );
1425
1426 // Update solution ( x += U_*C_'*r_ )
1427 RCP<MV> update = MVT::Clone( *problem_->getCurrLHSVec(), 1 );
1428 MVT::MvInit( *update, 0.0 );
1429 MVT::MvTimesMatAddMv( one, *Utmp, Ctr, one, *update );
1430 problem_->updateSolution( update, true );
1431
1432 // Update residual norm ( r -= C_*C_'*r_ )
1433 MVT::MvTimesMatAddMv( -one, *Ctmp, Ctr, one, *r_ );
1434
1435 // We recycled space from previous call
1436 prime_iterations = 0;
1437
1438 }
1439 else {
1440
1441 // Do one cycle of Gmres to "prime the pump" if there is no subspace to recycle
1442 printer_->stream(Debug) << " No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
1443
1444 Teuchos::ParameterList primeList;
1445
1446 // Tell the block solver that the block size is one.
1447 primeList.set("Num Blocks",numBlocks_);
1448 primeList.set("Recycled Blocks",0);
1449
1450 // Create GCRODR iterator object to perform one cycle of GMRES.
1451 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_prime_iter;
1452 gcrodr_prime_iter = rcp( new GCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,primeList) );
1453
1454 // Create the first block in the current Krylov basis (residual).
1455 problem_->computeCurrPrecResVec( &*r_ );
1456 index.resize( 1 ); index[0] = 0;
1457 RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1458 MVT::SetBlock(*r_,index,*v0); // V(:,0) = r
1459
1460 // Set the new state and initialize the solver.
1462 index.resize( numBlocks_+1 );
1463 for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1464 newstate.V = MVT::CloneViewNonConst( *V_, index );
1465 newstate.U = Teuchos::null;
1466 newstate.C = Teuchos::null;
1467 newstate.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, recycledBlocks_+1, recycledBlocks_+1 ) );
1468 newstate.B = Teuchos::null;
1469 newstate.curDim = 0;
1470 gcrodr_prime_iter->initialize(newstate);
1471
1472 // Perform one cycle of GMRES
1473 bool primeConverged = false;
1474 try {
1475 gcrodr_prime_iter->iterate();
1476
1477 // Check convergence first
1478 if ( convTest_->getStatus() == Passed ) {
1479 // we have convergence
1480 primeConverged = true;
1481 }
1482 }
1483 catch (const GCRODRIterOrthoFailure &e) {
1484 // Try to recover the most recent least-squares solution
1485 gcrodr_prime_iter->updateLSQR( gcrodr_prime_iter->getCurSubspaceDim() );
1486
1487 // Check to see if the most recent least-squares solution yielded convergence.
1488 sTest_->checkStatus( &*gcrodr_prime_iter );
1489 if (convTest_->getStatus() == Passed)
1490 primeConverged = true;
1491 }
1492 catch (const StatusTestNaNError& e) {
1493 // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1494 achievedTol_ = MT::one();
1495 Teuchos::RCP<MV> X = problem_->getLHS();
1496 MVT::MvInit( *X, SCT::zero() );
1497 printer_->stream(Warnings) << "Belos::GCRODRSolMgr::solve(): Warning! NaN has been detected!"
1498 << std::endl;
1499 return Unconverged;
1500 }
1501 catch (const std::exception &e) {
1502 printer_->stream(Errors) << "Error! Caught exception in GCRODRIter::iterate() at iteration "
1503 << gcrodr_prime_iter->getNumIters() << std::endl
1504 << e.what() << std::endl;
1505 throw;
1506 }
1507 // Record number of iterations in generating initial recycle spacec
1508 prime_iterations = gcrodr_prime_iter->getNumIters();
1509
1510 // Update the linear problem.
1511 RCP<MV> update = gcrodr_prime_iter->getCurrentUpdate();
1512 problem_->updateSolution( update, true );
1513
1514 // Get the state.
1515 newstate = gcrodr_prime_iter->getState();
1516 int p = newstate.curDim;
1517
1518 // Compute harmonic Ritz vectors
1519 // NOTE: The storage for the harmonic Ritz vectors (PP) is made one column larger
1520 // just in case we split a complex conjugate pair.
1521 // NOTE: Generate a recycled subspace only if we have enough vectors. If we converged
1522 // too early, move on to the next linear system and try to generate a subspace again.
1523 if (recycledBlocks_ < p+1) {
1524 int info = 0;
1525 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PPtmp = rcp (new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *PP_, p, recycledBlocks_+1 ) );
1526 // getHarmonicVecs1 assumes PP has recycledBlocks_+1 columns available
1527 keff = getHarmonicVecs1( p, *newstate.H, *PPtmp );
1528 // Hereafter, only keff columns of PP are needed
1529 PPtmp = rcp (new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *PP_, p, keff ) );
1530 // Now get views into C, U, V
1531 index.resize(keff);
1532 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1533 RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1534 RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1535 RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1536 index.resize(p);
1537 for (int ii=0; ii < p; ++ii) { index[ii] = ii; }
1538 RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1539
1540 // Form U (the subspace to recycle)
1541 // U = newstate.V(:,1:p) * PP;
1542 MVT::MvTimesMatAddMv( one, *Vtmp, *PPtmp, zero, *U1tmp );
1543
1544 // Form orthonormalized C and adjust U so that C = A*U
1545
1546 // First, compute [Q, R] = qr(H*P);
1547
1548 // Step #1: Form HP = H*P
1549 Teuchos::SerialDenseMatrix<int,ScalarType> Htmp( Teuchos::View, *H2_, p+1, p, recycledBlocks_+1,recycledBlocks_+1);
1550 Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *HP_, p+1, keff );
1551 HPtmp.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, Htmp, *PPtmp, zero );
1552
1553 // Step #1.5: Perform workspace size query for QR
1554 // factorization of HP. On input, lwork must be -1.
1555 // _GEQRF will put the workspace size in work_[0].
1556 int lwork = -1;
1557 tau_.resize (keff);
1558 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1559 HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1560 TEUCHOS_TEST_FOR_EXCEPTION(
1561 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve:"
1562 " LAPACK's _GEQRF failed to compute a workspace size.");
1563
1564 // Step #2: Compute QR factorization of HP
1565 //
1566 // NOTE (mfh 17 Apr 2014) LAPACK promises that the value of
1567 // work_[0] after the workspace query will fit in int. This
1568 // justifies the cast. We call real() first because
1569 // static_cast from std::complex to int doesn't work.
1570 lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0])));
1571 work_.resize (lwork); // Allocate workspace for the QR factorization
1572 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1573 HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1574 TEUCHOS_TEST_FOR_EXCEPTION(
1575 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve:"
1576 " LAPACK's _GEQRF failed to compute a QR factorization.");
1577
1578 // Step #3: Explicitly construct Q and R factors
1579 // NOTE: The upper triangular part of HP is copied into R and HP becomes Q.
1580 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff, keff );
1581 for (int ii = 0; ii < keff; ++ii) {
1582 for (int jj = ii; jj < keff; ++jj) {
1583 Rtmp(ii,jj) = HPtmp(ii,jj);
1584 }
1585 }
1586 // NOTE (mfh 17 Apr 2014): Teuchos::LAPACK's wrapper for
1587 // UNGQR dispatches to the correct Scalar-specific routine.
1588 // It calls {S,D}ORGQR if Scalar is real, and {C,Z}UNGQR if
1589 // Scalar is complex.
1590 lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
1591 HPtmp.values (), HPtmp.stride (), &tau_[0], &work_[0],
1592 lwork, &info);
1593 TEUCHOS_TEST_FOR_EXCEPTION(
1594 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1595 "LAPACK's _UNGQR failed to construct the Q factor.");
1596
1597 // Now we have [Q,R] = qr(H*P)
1598
1599 // Now compute C = V(:,1:p+1) * Q
1600 index.resize (p + 1);
1601 for (int ii = 0; ii < (p+1); ++ii) {
1602 index[ii] = ii;
1603 }
1604 Vtmp = MVT::CloneView( *V_, index ); // need new view into V (p+1 vectors now; needed p above)
1605 MVT::MvTimesMatAddMv( one, *Vtmp, HPtmp, zero, *Ctmp );
1606
1607 // Finally, compute U = U*R^{-1}.
1608 // This unfortuntely requires me to form R^{-1} explicitly and execute U = U * R^{-1}, as
1609 // backsolve capabilities don't exist in the Belos::MultiVec class
1610
1611 // Step #1: First, compute LU factorization of R
1612 ipiv_.resize(Rtmp.numRows());
1613 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1614 TEUCHOS_TEST_FOR_EXCEPTION(
1615 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1616 "LAPACK's _GETRF failed to compute an LU factorization.");
1617
1618 // FIXME (mfh 17 Apr 2014) We have to compute the explicit
1619 // inverse of R here because Belos::MultiVecTraits doesn't
1620 // have a triangular solve (where the triangular matrix is
1621 // globally replicated and the "right-hand side" is the
1622 // distributed MultiVector).
1623
1624 // Step #2: Form inv(R)
1625 lwork = Rtmp.numRows();
1626 work_.resize(lwork);
1627 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1628 TEUCHOS_TEST_FOR_EXCEPTION(
1629 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1630 "LAPACK's _GETRI failed to invert triangular matrix.");
1631
1632 // Step #3: Let U = U * R^{-1}
1633 MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
1634
1635 printer_->stream(Debug)
1636 << " Generated recycled subspace using RHS index " << currIdx[0]
1637 << " of dimension " << keff << std::endl << std::endl;
1638
1639 } // if (recycledBlocks_ < p+1)
1640
1641 // Return to outer loop if the priming solve converged, set the next linear system.
1642 if (primeConverged) {
1643 // Inform the linear problem that we are finished with this block linear system.
1644 problem_->setCurrLS();
1645
1646 // Update indices for the linear systems to be solved.
1647 numRHS2Solve -= 1;
1648 if (numRHS2Solve > 0) {
1649 currIdx[0]++;
1650 problem_->setLSIndex (currIdx); // Set the next indices
1651 }
1652 else {
1653 currIdx.resize (numRHS2Solve);
1654 }
1655
1656 continue;
1657 }
1658 } // if (keff > 0) ...
1659
1660 // Prepare for the Gmres iterations with the recycled subspace.
1661
1662 // Set the current number of recycled blocks and subspace dimension with the GCRO-DR iteration.
1663 gcrodr_iter->setSize( keff, numBlocks_ );
1664
1665 // Reset the number of iterations.
1666 gcrodr_iter->resetNumIters(prime_iterations);
1667
1668 // Reset the number of calls that the status test output knows about.
1669 outputTest_->resetNumCalls();
1670
1671 // Compute the residual after the priming solve, it will be the first block in the current Krylov basis.
1672 problem_->computeCurrPrecResVec( &*r_ );
1673 index.resize( 1 ); index[0] = 0;
1674 RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1675 MVT::SetBlock(*r_,index,*v0); // V(:,0) = r
1676
1677 // Set the new state and initialize the solver.
1679 index.resize( numBlocks_+1 );
1680 for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1681 newstate.V = MVT::CloneViewNonConst( *V_, index );
1682 index.resize( keff );
1683 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1684 newstate.C = MVT::CloneViewNonConst( *C_, index );
1685 newstate.U = MVT::CloneViewNonConst( *U_, index );
1686 newstate.B = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, keff, numBlocks_, 0, keff ) );
1687 newstate.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, keff, keff ) );
1688 newstate.curDim = 0;
1689 gcrodr_iter->initialize(newstate);
1690
1691 // variables needed for inner loop
1692 int numRestarts = 0;
1693 while(1) {
1694
1695 // tell gcrodr_iter to iterate
1696 try {
1697 gcrodr_iter->iterate();
1698
1700 //
1701 // check convergence first
1702 //
1704 if ( convTest_->getStatus() == Passed ) {
1705 // we have convergence
1706 break; // break from while(1){gcrodr_iter->iterate()}
1707 }
1709 //
1710 // check for maximum iterations
1711 //
1713 else if ( maxIterTest_->getStatus() == Passed ) {
1714 // we don't have convergence
1715 isConverged = false;
1716 break; // break from while(1){gcrodr_iter->iterate()}
1717 }
1719 //
1720 // check for restarting, i.e. the subspace is full
1721 //
1723 else if ( gcrodr_iter->getCurSubspaceDim() == gcrodr_iter->getMaxSubspaceDim() ) {
1724
1725 // Update the recycled subspace even if we have hit the maximum number of restarts.
1726
1727 // Update the linear problem.
1728 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1729 problem_->updateSolution( update, true );
1730
1731 buildRecycleSpace2(gcrodr_iter);
1732
1733 printer_->stream(Debug)
1734 << " Generated new recycled subspace using RHS index "
1735 << currIdx[0] << " of dimension " << keff << std::endl
1736 << std::endl;
1737
1738 // NOTE: If we have hit the maximum number of restarts then quit
1739 if (numRestarts >= maxRestarts_) {
1740 isConverged = false;
1741 break; // break from while(1){gcrodr_iter->iterate()}
1742 }
1743 numRestarts++;
1744
1745 printer_->stream(Debug)
1746 << " Performing restart number " << numRestarts << " of "
1747 << maxRestarts_ << std::endl << std::endl;
1748
1749 // Create the restart vector (first block in the current Krylov basis)
1750 problem_->computeCurrPrecResVec( &*r_ );
1751 index.resize( 1 ); index[0] = 0;
1752 RCP<MV> v00 = MVT::CloneViewNonConst( *V_, index );
1753 MVT::SetBlock(*r_,index,*v00); // V(:,0) = r
1754
1755 // Set the new state and initialize the solver.
1756 GCRODRIterState<ScalarType,MV> restartState;
1757 index.resize( numBlocks_+1 );
1758 for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1759 restartState.V = MVT::CloneViewNonConst( *V_, index );
1760 index.resize( keff );
1761 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1762 restartState.U = MVT::CloneViewNonConst( *U_, index );
1763 restartState.C = MVT::CloneViewNonConst( *C_, index );
1764 restartState.B = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, keff, numBlocks_, 0, keff ) );
1765 restartState.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, keff, keff ) );
1766 restartState.curDim = 0;
1767 gcrodr_iter->initialize(restartState);
1768
1769
1770 } // end of restarting
1771
1773 //
1774 // we returned from iterate(), but none of our status tests Passed.
1775 // something is wrong, and it is probably our fault.
1776 //
1778
1779 else {
1780 TEUCHOS_TEST_FOR_EXCEPTION(
1781 true, std::logic_error, "Belos::GCRODRSolMgr::solve: "
1782 "Invalid return from GCRODRIter::iterate().");
1783 }
1784 }
1785 catch (const GCRODRIterOrthoFailure &e) {
1786 // Try to recover the most recent least-squares solution
1787 gcrodr_iter->updateLSQR( gcrodr_iter->getCurSubspaceDim() );
1788
1789 // Check to see if the most recent least-squares solution yielded convergence.
1790 sTest_->checkStatus( &*gcrodr_iter );
1791 if (convTest_->getStatus() != Passed)
1792 isConverged = false;
1793 break;
1794 }
1795 catch (const std::exception& e) {
1796 printer_->stream(Errors)
1797 << "Error! Caught exception in GCRODRIter::iterate() at iteration "
1798 << gcrodr_iter->getNumIters() << std::endl << e.what() << std::endl;
1799 throw;
1800 }
1801 }
1802
1803 // Compute the current solution.
1804 // Update the linear problem.
1805 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1806 problem_->updateSolution( update, true );
1807
1808 // Inform the linear problem that we are finished with this block linear system.
1809 problem_->setCurrLS();
1810
1811 // If we didn't build a recycle space this solve but ran at least k iterations,
1812 // force build of new recycle space
1813
1814 if (!builtRecycleSpace_) {
1815 buildRecycleSpace2(gcrodr_iter);
1816 printer_->stream(Debug)
1817 << " Generated new recycled subspace using RHS index " << currIdx[0]
1818 << " of dimension " << keff << std::endl << std::endl;
1819 }
1820
1821 // Update indices for the linear systems to be solved.
1822 numRHS2Solve -= 1;
1823 if (numRHS2Solve > 0) {
1824 currIdx[0]++;
1825 problem_->setLSIndex (currIdx); // Set the next indices
1826 }
1827 else {
1828 currIdx.resize (numRHS2Solve);
1829 }
1830 } // while (numRHS2Solve > 0)
1831 }
1832
1833 sTest_->print (printer_->stream (FinalSummary)); // print final summary
1834
1835 // print timing information
1836#ifdef BELOS_TEUCHOS_TIME_MONITOR
1837 // Calling summarize() can be expensive, so don't call unless the
1838 // user wants to print out timing details. summarize() will do all
1839 // the work even if it's passed a "black hole" output stream.
1840 if (verbosity_ & TimingDetails)
1841 Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1842#endif // BELOS_TEUCHOS_TIME_MONITOR
1843
1844 // get iteration information for this solve
1845 numIters_ = maxIterTest_->getNumIters ();
1846
1847 // Save the convergence test value ("achieved tolerance") for this
1848 // solve. This solver (unlike BlockGmresSolMgr) always has two
1849 // residual norm status tests: an explicit and an implicit test.
1850 // The master convergence test convTest_ is a SEQ combo of the
1851 // implicit resp. explicit tests. If the implicit test never
1852 // passes, then the explicit test won't ever be executed. This
1853 // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1854 // with this case by using the values returned by
1855 // impConvTest_->getTestValue().
1856 {
1857 const std::vector<MagnitudeType>* pTestValues = expConvTest_->getTestValue();
1858 if (pTestValues == NULL || pTestValues->size() < 1) {
1859 pTestValues = impConvTest_->getTestValue();
1860 }
1861 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1862 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1863 "method returned NULL. Please report this bug to the Belos developers.");
1864 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1865 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1866 "method returned a vector of length zero. Please report this bug to the "
1867 "Belos developers.");
1868
1869 // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1870 // achieved tolerances for all vectors in the current solve(), or
1871 // just for the vectors from the last deflation?
1872 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1873 }
1874
1875 return isConverged ? Converged : Unconverged; // return from solve()
1876}
1877
1878// Given existing recycle space and Krylov space, build new recycle space
1879template<class ScalarType, class MV, class OP>
1881
1882 MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
1883 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1884
1885 std::vector<MagnitudeType> d(keff);
1886 std::vector<ScalarType> dscalar(keff);
1887 std::vector<int> index(numBlocks_+1);
1888
1889 // Get the state
1890 GCRODRIterState<ScalarType,MV> oldState = gcrodr_iter->getState();
1891 int p = oldState.curDim;
1892
1893 // insufficient new information to update recycle space
1894 if (p<1) return;
1895
1896 // Take the norm of the recycled vectors
1897 {
1898 index.resize(keff);
1899 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1900 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1901 d.resize(keff);
1902 dscalar.resize(keff);
1903 MVT::MvNorm( *Utmp, d );
1904 for (int i=0; i<keff; ++i) {
1905 d[i] = one / d[i];
1906 dscalar[i] = (ScalarType)d[i];
1907 }
1908 MVT::MvScale( *Utmp, dscalar );
1909 }
1910
1911 // Get view into current "full" upper Hessnburg matrix
1912 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H2tmp = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, p+keff+1, p+keff ) );
1913
1914 // Insert D into the leading keff x keff block of H2
1915 for (int i=0; i<keff; ++i) {
1916 (*H2tmp)(i,i) = d[i];
1917 }
1918
1919 // Compute the harmoic Ritz pairs for the generalized eigenproblem
1920 // getHarmonicVecs2 assumes PP has recycledBlocks_+1 columns available
1921 int keff_new;
1922 {
1923 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p+keff, recycledBlocks_+1 );
1924 keff_new = getHarmonicVecs2( keff, p, *H2tmp, oldState.V, PPtmp );
1925 }
1926
1927 // Code to form new U, C
1928 // U = [U V(:,1:p)] * P; (in two steps)
1929
1930 // U(:,1:keff) = matmul(U(:,1:keff_old),PP(1:keff_old,1:keff)) (step 1)
1931 Teuchos::RCP<MV> U1tmp;
1932 {
1933 index.resize( keff );
1934 for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1935 Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1936 index.resize( keff_new );
1937 for (int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1938 U1tmp = MVT::CloneViewNonConst( *U1_, index );
1939 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, keff, keff_new );
1940 MVT::MvTimesMatAddMv( one, *Utmp, PPtmp, zero, *U1tmp );
1941 }
1942
1943 // U(:,1:keff) = U(:,1:keff) + matmul(V(:,1:m-k),PP(keff_old+1:m-k+keff_old,1:keff)) (step 2)
1944 {
1945 index.resize(p);
1946 for (int ii=0; ii < p; ii++) { index[ii] = ii; }
1947 Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1948 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p, keff_new, keff );
1949 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *U1tmp );
1950 }
1951
1952 // Form HP = H*P
1953 Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *HP_, p+keff+1, keff_new );
1954 {
1955 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p+keff, keff_new );
1956 HPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,*H2tmp,PPtmp,zero);
1957 }
1958
1959 // Workspace size query for QR factorization of HP (the worksize will be placed in work_[0])
1960 int info = 0, lwork = -1;
1961 tau_.resize (keff_new);
1962 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1963 HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1964 TEUCHOS_TEST_FOR_EXCEPTION(
1965 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1966 "LAPACK's _GEQRF failed to compute a workspace size.");
1967
1968 // NOTE (mfh 18 Apr 2014) LAPACK promises that the value of work_[0]
1969 // after the workspace query will fit in int. This justifies the
1970 // cast. We call real() first because static_cast from std::complex
1971 // to int doesn't work.
1972 lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0])));
1973 work_.resize (lwork); // Allocate workspace for the QR factorization
1974 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1975 HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1976 TEUCHOS_TEST_FOR_EXCEPTION(
1977 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1978 "LAPACK's _GEQRF failed to compute a QR factorization.");
1979
1980 // Explicitly construct Q and R factors
1981 // NOTE: The upper triangular part of HP is copied into R and HP becomes Q.
1982 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff_new, keff_new );
1983 for(int i=0;i<keff_new;i++) { for(int j=i;j<keff_new;j++) Rtmp(i,j) = HPtmp(i,j); }
1984
1985 // NOTE (mfh 18 Apr 2014): Teuchos::LAPACK's wrapper for UNGQR
1986 // dispatches to the correct Scalar-specific routine. It calls
1987 // {S,D}ORGQR if Scalar is real, and {C,Z}UNGQR if Scalar is
1988 // complex.
1989 lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
1990 HPtmp.values (), HPtmp.stride (), &tau_[0], &work_[0],
1991 lwork, &info);
1992 TEUCHOS_TEST_FOR_EXCEPTION(
1993 info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1994 "LAPACK's _UNGQR failed to construct the Q factor.");
1995
1996 // Form orthonormalized C and adjust U accordingly so that C = A*U
1997 // C = [C V] * Q;
1998
1999 // C(:,1:keff) = matmul(C(:,1:keff_old),QQ(1:keff_old,1:keff))
2000 {
2001 Teuchos::RCP<MV> C1tmp;
2002 {
2003 index.resize(keff);
2004 for (int i=0; i < keff; i++) { index[i] = i; }
2005 Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
2006 index.resize(keff_new);
2007 for (int i=0; i < keff_new; i++) { index[i] = i; }
2008 C1tmp = MVT::CloneViewNonConst( *C1_, index );
2009 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *HP_, keff, keff_new );
2010 MVT::MvTimesMatAddMv( one, *Ctmp, PPtmp, zero, *C1tmp );
2011 }
2012 // Now compute C += V(:,1:p+1) * Q
2013 {
2014 index.resize( p+1 );
2015 for (int i=0; i < p+1; ++i) { index[i] = i; }
2016 Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
2017 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *HP_, p+1, keff_new, keff, 0 );
2018 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *C1tmp );
2019 }
2020 }
2021
2022 // C_ = C1_; (via a swap)
2023 std::swap(C_, C1_);
2024
2025 // Finally, compute U_ = U_*R^{-1}
2026 // First, compute LU factorization of R
2027 ipiv_.resize(Rtmp.numRows());
2028 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
2029 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
2030
2031 // Now, form inv(R)
2032 lwork = Rtmp.numRows();
2033 work_.resize(lwork);
2034 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
2035 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
2036
2037 {
2038 index.resize(keff_new);
2039 for (int i=0; i < keff_new; i++) { index[i] = i; }
2040 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
2041 MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
2042 }
2043
2044 // Set the current number of recycled blocks and subspace dimension with the GCRO-DR iteration.
2045 if (keff != keff_new) {
2046 keff = keff_new;
2047 gcrodr_iter->setSize( keff, numBlocks_ );
2048 // Important to zero this out before next cyle
2049 Teuchos::SerialDenseMatrix<int,ScalarType> b1( Teuchos::View, *H2_, recycledBlocks_+2, 1, 0, recycledBlocks_ );
2050 b1.putScalar(zero);
2051 }
2052
2053}
2054
2055
2056// Compute the harmonic eigenpairs of the projected, dense system.
2057template<class ScalarType, class MV, class OP>
2059 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2060 Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2061 int i, j;
2062 bool xtraVec = false;
2063 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2064
2065 // Real and imaginary eigenvalue components
2066 std::vector<MagnitudeType> wr(m), wi(m);
2067
2068 // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
2069 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m,m,false);
2070
2071 // Magnitude of harmonic Ritz values
2072 std::vector<MagnitudeType> w(m);
2073
2074 // Sorted order of harmonic Ritz values, also used for GEEV
2075 std::vector<int> iperm(m);
2076
2077 // Output info
2078 int info = 0;
2079
2080 // Set flag indicating recycle space has been generated this solve
2081 builtRecycleSpace_ = true;
2082
2083 // Solve linear system: H_m^{-H}*e_m
2084 Teuchos::SerialDenseMatrix<int, ScalarType> HHt( HH, Teuchos::TRANS );
2085 Teuchos::SerialDenseVector<int, ScalarType> e_m( m );
2086 e_m[m-1] = one;
2087 lapack.GESV(m, 1, HHt.values(), HHt.stride(), &iperm[0], e_m.values(), e_m.stride(), &info);
2088 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
2089
2090 // Compute H_m + d*H_m^{-H}*e_m*e_m^H
2091 ScalarType d = HH(m, m-1) * HH(m, m-1);
2092 Teuchos::SerialDenseMatrix<int, ScalarType> harmHH( Teuchos::Copy, HH, m, m );
2093 for( i=0; i<m; ++i )
2094 harmHH(i, m-1) += d * e_m[i];
2095
2096 // Revise to do query for optimal workspace first
2097 // Create simple storage for the left eigenvectors, which we don't care about.
2098 const int ldvl = 1;
2099 ScalarType* vl = 0;
2100
2101 // Size of workspace and workspace for GEEV
2102 int lwork = -1;
2103 std::vector<ScalarType> work(1);
2104 std::vector<MagnitudeType> rwork(2*m);
2105
2106 // First query GEEV for the optimal workspace size
2107 lapack.GEEV('N', 'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2108 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
2109
2110 lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
2111 work.resize( lwork );
2112
2113 lapack.GEEV('N', 'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2114 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
2115 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
2116
2117 // Construct magnitude of each harmonic Ritz value
2118 for( i=0; i<m; ++i )
2119 w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
2120
2121 // Construct magnitude of each harmonic Ritz value
2122 this->sort(w, m, iperm);
2123
2124 const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2125
2126 // Select recycledBlocks_ smallest eigenvectors
2127 for( i=0; i<recycledBlocks_; ++i ) {
2128 for( j=0; j<m; j++ ) {
2129 PP(j,i) = vr(j,iperm[i]);
2130 }
2131 }
2132
2133 if(!scalarTypeIsComplex) {
2134
2135 // Determine exact size for PP (i.e., determine if we need to store an additional vector)
2136 if (wi[iperm[recycledBlocks_-1]] != 0.0) {
2137 int countImag = 0;
2138 for ( i=0; i<recycledBlocks_; ++i ) {
2139 if (wi[iperm[i]] != 0.0)
2140 countImag++;
2141 }
2142 // Check to see if this count is even or odd:
2143 if (countImag % 2)
2144 xtraVec = true;
2145 }
2146
2147 if (xtraVec) { // we need to store one more vector
2148 if (wi[iperm[recycledBlocks_-1]] > 0.0) { // I picked the "real" component
2149 for( j=0; j<m; ++j ) { // so get the "imag" component
2150 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
2151 }
2152 }
2153 else { // I picked the "imag" component
2154 for( j=0; j<m; ++j ) { // so get the "real" component
2155 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
2156 }
2157 }
2158 }
2159
2160 }
2161
2162 // Return whether we needed to store an additional vector
2163 if (xtraVec) {
2164 return recycledBlocks_+1;
2165 }
2166 else {
2167 return recycledBlocks_;
2168 }
2169
2170}
2171
2172// Compute the harmonic eigenpairs of the projected, dense system.
2173template<class ScalarType, class MV, class OP>
2175 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2176 const Teuchos::RCP<const MV>& VV,
2177 Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2178 int i, j;
2179 int m2 = HH.numCols();
2180 bool xtraVec = false;
2181 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2182 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
2183 std::vector<int> index;
2184
2185 // Real and imaginary eigenvalue components
2186 std::vector<MagnitudeType> wr(m2), wi(m2);
2187
2188 // Magnitude of harmonic Ritz values
2189 std::vector<MagnitudeType> w(m2);
2190
2191 // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
2192 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m2,m2,false);
2193
2194 // Sorted order of harmonic Ritz values
2195 std::vector<int> iperm(m2);
2196
2197 // Set flag indicating recycle space has been generated this solve
2198 builtRecycleSpace_ = true;
2199
2200 // Form matrices for generalized eigenproblem
2201
2202 // B = H2' * H2; Don't zero out matrix when constructing
2203 Teuchos::SerialDenseMatrix<int,ScalarType> B(m2,m2,false);
2204 B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,HH,HH,zero);
2205
2206 // A_tmp = | C'*U 0 |
2207 // | V_{m+1}'*U I |
2208 Teuchos::SerialDenseMatrix<int,ScalarType> A_tmp( keffloc+m+1, keffloc+m );
2209
2210 // A_tmp(1:keffloc,1:keffloc) = C' * U;
2211 index.resize(keffloc);
2212 for (i=0; i<keffloc; ++i) { index[i] = i; }
2213 Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
2214 Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
2215 Teuchos::SerialDenseMatrix<int,ScalarType> A11( Teuchos::View, A_tmp, keffloc, keffloc );
2216 MVT::MvTransMv( one, *Ctmp, *Utmp, A11 );
2217
2218 // A_tmp(keffloc+1:m-k+keffloc+1,1:keffloc) = V' * U;
2219 Teuchos::SerialDenseMatrix<int,ScalarType> A21( Teuchos::View, A_tmp, m+1, keffloc, keffloc );
2220 index.resize(m+1);
2221 for (i=0; i < m+1; i++) { index[i] = i; }
2222 Teuchos::RCP<const MV> Vp = MVT::CloneView( *VV, index );
2223 MVT::MvTransMv( one, *Vp, *Utmp, A21 );
2224
2225 // A_tmp(keffloc+1:m-k+keffloc,keffloc+1:m-k+keffloc) = eye(m-k);
2226 for( i=keffloc; i<keffloc+m; i++ ) {
2227 A_tmp(i,i) = one;
2228 }
2229
2230 // A = H2' * A_tmp;
2231 Teuchos::SerialDenseMatrix<int,ScalarType> A( m2, A_tmp.numCols() );
2232 A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, HH, A_tmp, zero );
2233
2234 // Compute k smallest harmonic Ritz pairs
2235 // SUBROUTINE DGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
2236 // ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO,
2237 // IHI, LSCALE, RSCALE, ABNRM, BBNRM, RCONDE,
2238 // RCONDV, WORK, LWORK, IWORK, BWORK, INFO )
2239 // MLP: 'SCALING' in DGGEVX generates incorrect eigenvalues. Therefore, only permuting
2240 char balanc='P', jobvl='N', jobvr='V', sense='N';
2241 int ld = A.numRows();
2242 int lwork = 6*ld;
2243 int ldvl = ld, ldvr = ld;
2244 int info = 0,ilo = 0,ihi = 0;
2245 MagnitudeType abnrm = 0.0, bbnrm = 0.0;
2246 ScalarType *vl = 0; // This is never referenced by dggevx if jobvl == 'N'
2247 std::vector<ScalarType> beta(ld);
2248 std::vector<ScalarType> work(lwork);
2249 std::vector<MagnitudeType> rwork(lwork);
2250 std::vector<MagnitudeType> lscale(ld), rscale(ld);
2251 std::vector<MagnitudeType> rconde(ld), rcondv(ld);
2252 std::vector<int> iwork(ld+6);
2253 int *bwork = 0; // If sense == 'N', bwork is never referenced
2254 //lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
2255 // &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2256 // &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &iwork[0], bwork, &info);
2257 lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
2258 &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2259 &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
2260 &iwork[0], bwork, &info);
2261 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
2262
2263 // Construct magnitude of each harmonic Ritz value
2264 // NOTE : Forming alpha/beta *should* be okay here, given assumptions on construction of matrix pencil above
2265 for( i=0; i<ld; i++ ) {
2266 w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot (wr[i]*wr[i] + wi[i]*wi[i]) /
2267 Teuchos::ScalarTraits<ScalarType>::magnitude (beta[i]);
2268 }
2269
2270 // Construct magnitude of each harmonic Ritz value
2271 this->sort(w,ld,iperm);
2272
2273 const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2274
2275 // Select recycledBlocks_ smallest eigenvectors
2276 for( i=0; i<recycledBlocks_; i++ ) {
2277 for( j=0; j<ld; j++ ) {
2278 PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
2279 }
2280 }
2281
2282 if(!scalarTypeIsComplex) {
2283
2284 // Determine exact size for PP (i.e., determine if we need to store an additional vector)
2285 if (wi[iperm[ld-recycledBlocks_]] != 0.0) {
2286 int countImag = 0;
2287 for ( i=ld-recycledBlocks_; i<ld; i++ ) {
2288 if (wi[iperm[i]] != 0.0)
2289 countImag++;
2290 }
2291 // Check to see if this count is even or odd:
2292 if (countImag % 2)
2293 xtraVec = true;
2294 }
2295
2296 if (xtraVec) { // we need to store one more vector
2297 if (wi[iperm[ld-recycledBlocks_]] > 0.0) { // I picked the "real" component
2298 for( j=0; j<ld; j++ ) { // so get the "imag" component
2299 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
2300 }
2301 }
2302 else { // I picked the "imag" component
2303 for( j=0; j<ld; j++ ) { // so get the "real" component
2304 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
2305 }
2306 }
2307 }
2308
2309 }
2310
2311 // Return whether we needed to store an additional vector
2312 if (xtraVec) {
2313 return recycledBlocks_+1;
2314 }
2315 else {
2316 return recycledBlocks_;
2317 }
2318
2319}
2320
2321
2322// This method sorts list of n floating-point numbers and return permutation vector
2323template<class ScalarType, class MV, class OP>
2324void GCRODRSolMgr<ScalarType,MV,OP,true>::sort(std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm) {
2325 int l, r, j, i, flag;
2326 int RR2;
2327 MagnitudeType dRR, dK;
2328
2329 // Initialize the permutation vector.
2330 for(j=0;j<n;j++)
2331 iperm[j] = j;
2332
2333 if (n <= 1) return;
2334
2335 l = n / 2 + 1;
2336 r = n - 1;
2337 l = l - 1;
2338 dRR = dlist[l - 1];
2339 dK = dlist[l - 1];
2340
2341 RR2 = iperm[l - 1];
2342 while (r != 0) {
2343 j = l;
2344 flag = 1;
2345
2346 while (flag == 1) {
2347 i = j;
2348 j = j + j;
2349
2350 if (j > r + 1)
2351 flag = 0;
2352 else {
2353 if (j < r + 1)
2354 if (dlist[j] > dlist[j - 1]) j = j + 1;
2355
2356 if (dlist[j - 1] > dK) {
2357 dlist[i - 1] = dlist[j - 1];
2358 iperm[i - 1] = iperm[j - 1];
2359 }
2360 else {
2361 flag = 0;
2362 }
2363 }
2364 }
2365 dlist[i - 1] = dRR;
2366 iperm[i - 1] = RR2;
2367
2368 if (l == 1) {
2369 dRR = dlist [r];
2370 RR2 = iperm[r];
2371 dK = dlist[r];
2372 dlist[r] = dlist[0];
2373 iperm[r] = iperm[0];
2374 r = r - 1;
2375 }
2376 else {
2377 l = l - 1;
2378 dRR = dlist[l - 1];
2379 RR2 = iperm[l - 1];
2380 dK = dlist[l - 1];
2381 }
2382 }
2383 dlist[0] = dRR;
2384 iperm[0] = RR2;
2385}
2386
2387
2388template<class ScalarType, class MV, class OP>
2390 std::ostringstream out;
2391 out << "Belos::GCRODRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
2392 out << "{";
2393 out << "Ortho Type: \"" << orthoType_ << "\"";
2394 out << ", Num Blocks: " <<numBlocks_;
2395 out << ", Num Recycle Blocks: " << recycledBlocks_;
2396 out << ", Max Restarts: " << maxRestarts_;
2397 out << "}";
2398 return out.str ();
2399}
2400
2401} // namespace Belos
2402
2403#endif /* BELOS_GCRODR_SOLMGR_HPP */
Belos concrete class for performing the block, flexible GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Belos concrete class for performing the GCRO-DR iteration.
Class which describes the linear problem to be solved by the iterative solver.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class which describes the basic interface for a solver manager.
Belos::StatusTest for logically combining several status tests.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
Collection of types and exceptions used within the Belos solvers.
BelosError(const std::string &what_arg)
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
void setDepTol(const MagnitudeType dep_tol)
Set parameter for re-orthogonalization threshhold.
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
This class implements the GCRODR iteration, where a single-stdvector Krylov subspace is constructed....
GCRODRSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call.
GCRODRSolMgrLAPACKFailure(const std::string &what_arg)
GCRODRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i....
GCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
GCRODRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorm...
GCRODRSolMgrOrthoFailure(const std::string &what_arg)
GCRODRSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspa...
GCRODRSolMgrRecyclingFailure(const std::string &what_arg)
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
GCRODRSolMgr()
Empty constructor for GCRODRSolMgr. This constructor takes no arguments and sets the default values f...
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
Implementation of the GCRODR (Recycling GMRES) iterative linear solver.
GCRODRSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Traits class which defines basic operations on multivectors.
static void MvTransMv(const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Class which defines basic traits for the operator type.
Enumeration of all valid Belos (Mat)OrthoManager classes.
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
Return an instance of the specified MatOrthoManager subclass.
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
A class for extending the status testing capabilities of Belos via logical combinations.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
An implementation of StatusTestResNorm using a family of residual norms.
A Belos::StatusTest class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
@ StatusTestDetails
@ FinalSummary
@ TimingDetails
ReturnType
Whether the Belos solve converged for all linear systems.
@ Unconverged
ScaleType
The type of scaling to use on the residual norm value.
ResetType
How to reset the solver.
@ RecycleSubspace
static const double convTol
Default convergence tolerance.
Structure to contain pointers to GCRODRIter state variables.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Teuchos::RCP< MV > V
The current Krylov basis.
int curDim
The current dimension of the reduction.

Generated on for Belos by doxygen 1.15.0