This file is indexed.

/usr/include/trilinos/BelosBlockCGSolMgr.hpp is in libtrilinos-belos-dev 12.10.1-3.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

   1
   2
   3
   4
   5
   6
   7
   8
   9
  10
  11
  12
  13
  14
  15
  16
  17
  18
  19
  20
  21
  22
  23
  24
  25
  26
  27
  28
  29
  30
  31
  32
  33
  34
  35
  36
  37
  38
  39
  40
  41
  42
  43
  44
  45
  46
  47
  48
  49
  50
  51
  52
  53
  54
  55
  56
  57
  58
  59
  60
  61
  62
  63
  64
  65
  66
  67
  68
  69
  70
  71
  72
  73
  74
  75
  76
  77
  78
  79
  80
  81
  82
  83
  84
  85
  86
  87
  88
  89
  90
  91
  92
  93
  94
  95
  96
  97
  98
  99
 100
 101
 102
 103
 104
 105
 106
 107
 108
 109
 110
 111
 112
 113
 114
 115
 116
 117
 118
 119
 120
 121
 122
 123
 124
 125
 126
 127
 128
 129
 130
 131
 132
 133
 134
 135
 136
 137
 138
 139
 140
 141
 142
 143
 144
 145
 146
 147
 148
 149
 150
 151
 152
 153
 154
 155
 156
 157
 158
 159
 160
 161
 162
 163
 164
 165
 166
 167
 168
 169
 170
 171
 172
 173
 174
 175
 176
 177
 178
 179
 180
 181
 182
 183
 184
 185
 186
 187
 188
 189
 190
 191
 192
 193
 194
 195
 196
 197
 198
 199
 200
 201
 202
 203
 204
 205
 206
 207
 208
 209
 210
 211
 212
 213
 214
 215
 216
 217
 218
 219
 220
 221
 222
 223
 224
 225
 226
 227
 228
 229
 230
 231
 232
 233
 234
 235
 236
 237
 238
 239
 240
 241
 242
 243
 244
 245
 246
 247
 248
 249
 250
 251
 252
 253
 254
 255
 256
 257
 258
 259
 260
 261
 262
 263
 264
 265
 266
 267
 268
 269
 270
 271
 272
 273
 274
 275
 276
 277
 278
 279
 280
 281
 282
 283
 284
 285
 286
 287
 288
 289
 290
 291
 292
 293
 294
 295
 296
 297
 298
 299
 300
 301
 302
 303
 304
 305
 306
 307
 308
 309
 310
 311
 312
 313
 314
 315
 316
 317
 318
 319
 320
 321
 322
 323
 324
 325
 326
 327
 328
 329
 330
 331
 332
 333
 334
 335
 336
 337
 338
 339
 340
 341
 342
 343
 344
 345
 346
 347
 348
 349
 350
 351
 352
 353
 354
 355
 356
 357
 358
 359
 360
 361
 362
 363
 364
 365
 366
 367
 368
 369
 370
 371
 372
 373
 374
 375
 376
 377
 378
 379
 380
 381
 382
 383
 384
 385
 386
 387
 388
 389
 390
 391
 392
 393
 394
 395
 396
 397
 398
 399
 400
 401
 402
 403
 404
 405
 406
 407
 408
 409
 410
 411
 412
 413
 414
 415
 416
 417
 418
 419
 420
 421
 422
 423
 424
 425
 426
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
//@HEADER
// ************************************************************************
//
//                 Belos: Block Linear Solvers Package
//                  Copyright 2004 Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ************************************************************************
//@HEADER

#ifndef BELOS_BLOCK_CG_SOLMGR_HPP
#define BELOS_BLOCK_CG_SOLMGR_HPP

/*! \file BelosBlockCGSolMgr.hpp
 *  \brief The Belos::BlockCGSolMgr provides a solver manager for the BlockCG linear solver.
*/

#include "BelosConfigDefs.hpp"
#include "BelosTypes.hpp"

#include "BelosLinearProblem.hpp"
#include "BelosSolverManager.hpp"

#include "BelosCGIter.hpp"
#include "BelosBlockCGIter.hpp"
#include "BelosDGKSOrthoManager.hpp"
#include "BelosICGSOrthoManager.hpp"
#include "BelosIMGSOrthoManager.hpp"
#include "BelosStatusTestMaxIters.hpp"
#include "BelosStatusTestGenResNorm.hpp"
#include "BelosStatusTestCombo.hpp"
#include "BelosStatusTestOutputFactory.hpp"
#include "BelosOutputManager.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_LAPACK.hpp"
#ifdef BELOS_TEUCHOS_TIME_MONITOR
#  include "Teuchos_TimeMonitor.hpp"
#endif
#if defined(HAVE_TEUCHOSCORE_CXX11)
#  include <type_traits>
#endif // defined(HAVE_TEUCHOSCORE_CXX11)
#include <algorithm>

/** \example BlockCG/BlockCGEpetraExFile.cpp
    This is an example of how to use the Belos::BlockCGSolMgr solver manager.
*/
/** \example BlockCG/BlockPrecCGEpetraExFile.cpp
    This is an example of how to use the Belos::BlockCGSolMgr solver manager with an Ifpack preconditioner.
*/

/*! \class Belos::BlockCGSolMgr
 *
 *  \brief The Belos::BlockCGSolMgr provides a powerful and fully-featured solver manager over the CG and BlockCG linear solver.

 \ingroup belos_solver_framework

 \author Heidi Thornquist, Chris Baker, and Teri Barth
 */

namespace Belos {

  //! @name BlockCGSolMgr Exceptions
  //@{

  /** \brief BlockCGSolMgrLinearProblemFailure is thrown when the linear problem is
   * not setup (i.e. setProblem() was not called) when solve() is called.
   *
   * This std::exception is thrown from the BlockCGSolMgr::solve() method.
   *
   */
  class BlockCGSolMgrLinearProblemFailure : public BelosError {public:
    BlockCGSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
    {}};

  /** \brief BlockCGSolMgrOrthoFailure is thrown when the orthogonalization manager is
   * unable to generate orthonormal columns from the initial basis vectors.
   *
   * This std::exception is thrown from the BlockCGSolMgr::solve() method.
   *
   */
  class BlockCGSolMgrOrthoFailure : public BelosError {public:
    BlockCGSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
    {}};


  template<class ScalarType, class MV, class OP,
           const bool lapackSupportsScalarType =
           Belos::Details::LapackSupportsScalar<ScalarType>::value>
  class BlockCGSolMgr :
    public Details::SolverManagerRequiresLapack<ScalarType,MV,OP>
  {
    static const bool requiresLapack =
      Belos::Details::LapackSupportsScalar<ScalarType>::value;
    typedef Details::SolverManagerRequiresLapack<ScalarType, MV, OP,
                                                 requiresLapack> base_type;

  public:
    BlockCGSolMgr () :
      base_type ()
    {}
    BlockCGSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >& problem,
                  const Teuchos::RCP<Teuchos::ParameterList>& pl) :
      base_type ()
    {}
    virtual ~BlockCGSolMgr () {}
  };


  // Partial specialization for ScalarType types for which
  // Teuchos::LAPACK has a valid implementation.  This contains the
  // actual working implementation of BlockCGSolMgr.
  template<class ScalarType, class MV, class OP>
  class BlockCGSolMgr<ScalarType, MV, OP, true> :
    public Details::SolverManagerRequiresLapack<ScalarType, MV, OP, true>
  {
    // This partial specialization is already chosen for those scalar types
    // that support lapack, so we don't need to have an additional compile-time
    // check that the scalar type is float/double/complex.
// #if defined(HAVE_TEUCHOSCORE_CXX11)
// #  if defined(HAVE_TEUCHOS_COMPLEX)
//     static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
//                    std::is_same<ScalarType, std::complex<double> >::value ||
//                    std::is_same<ScalarType, float>::value ||
//                    std::is_same<ScalarType, double>::value,
//                    "Belos::GCRODRSolMgr: ScalarType must be one of the four "
//                    "types (S,D,C,Z) supported by LAPACK.");
// #  else
//     static_assert (std::is_same<ScalarType, float>::value ||
//                    std::is_same<ScalarType, double>::value,
//                    "Belos::GCRODRSolMgr: ScalarType must be float or double.  "
//                    "Complex arithmetic support is currently disabled.  To "
//                    "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
// #  endif // defined(HAVE_TEUCHOS_COMPLEX)
// #endif // defined(HAVE_TEUCHOSCORE_CXX11)

  private:
    typedef MultiVecTraits<ScalarType,MV> MVT;
    typedef OperatorTraits<ScalarType,MV,OP> OPT;
    typedef Teuchos::ScalarTraits<ScalarType> SCT;
    typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
    typedef Teuchos::ScalarTraits<MagnitudeType> MT;

  public:

    //! @name Constructors/Destructor
    //@{

    /*! \brief Empty constructor for BlockCGSolMgr.
     * This constructor takes no arguments and sets the default values for the solver.
     * The linear problem must be passed in using setProblem() before solve() is called on this object.
     * The solver values can be changed using setParameters().
     */
     BlockCGSolMgr();

    /*! \brief Basic constructor for BlockCGSolMgr.
     *
     * This constructor accepts the LinearProblem to be solved in addition
     * to a parameter list of options for the solver manager. These options include the following:
     *   - "Block Size" - an \c int specifying the block size to be used by the underlying block
     *                    conjugate-gradient solver. Default: 1
     *   - "Adaptive Block Size" - a \c bool specifying whether the block size can be modified
     *                             throughout the solve. Default: true
     *   - "Maximum Iterations" - an \c int specifying the maximum number of iterations the
     *                            underlying solver is allowed to perform. Default: 1000
     *   - "Convergence Tolerance" - a \c MagnitudeType specifying the level that residual norms
     *                               must reach to decide convergence. Default: 1e-8.
     *   - "Orthogonalization" - a \c std::string specifying the desired orthogonalization:
     *                           DGKS ,ICGS, and IMGS. Default: "DGKS"
     *   - "Orthogonalization Constant" - a \c MagnitudeType used by DGKS orthogonalization to
     *                                    determine whether another step of classical Gram-Schmidt
     *                                    is necessary.  Default: -1 (use DGKS default)
     *   - "Verbosity" - a sum of MsgType specifying the verbosity. Default: Belos::Errors
     *   - "Output Style" - a OutputType specifying the style of output. Default: Belos::General
     *   - "Output Stream" - a reference-counted pointer to the output stream where all
     *                       solver output is sent.  Default: Teuchos::rcp(&std::cout,false)
     *   - "Output Frequency" - an \c int specifying how often convergence information should be
     *                          outputted.  Default: -1 (never)
     *   - "Show Maximum Residual Norm Only" - a \c bool specifying whether that only the maximum
     *                                         relative residual norm is printed if convergence
     *                                         information is printed. Default: false
     *   - "Timer Label" - a \c std::string to use as a prefix for the timer labels.  Default: "Belos"
     *		\param pl [in] ParameterList with construction information
     *			\htmlonly
     *			<iframe src="belos_BlockCG.xml" width=100% scrolling="no" frameborder="0">
     *			</iframe>
     *			<hr />
     *			\endhtmlonly
     */
    BlockCGSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
                   const Teuchos::RCP<Teuchos::ParameterList> &pl );

    //! Destructor.
    virtual ~BlockCGSolMgr() {};
    //@}

    //! @name Accessor methods
    //@{

    const LinearProblem<ScalarType,MV,OP>& getProblem() const {
      return *problem_;
    }

    /*! \brief Get a parameter list containing the valid parameters for this object.
     */
    Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;

    /*! \brief Get a parameter list containing the current parameters for this object.
     */
    Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }

    /*! \brief Return the timers for this object.
     *
     * The timers are ordered as follows:
     *   - time spent in solve() routine
     */
    Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
      return Teuchos::tuple(timerSolve_);
    }

    /// \brief Tolerance achieved by the last \c solve() invocation.
    ///
    /// This is the maximum over all right-hand sides' achieved
    /// convergence tolerances, and is set whether or not the solve
    /// actually managed to achieve the desired convergence tolerance.
    MagnitudeType achievedTol() const {
      return achievedTol_;
    }

    //! Get the iteration count for the most recent call to \c solve().
    int getNumIters() const {
      return numIters_;
    }

    /*! \brief Return whether a loss of accuracy was detected by this solver during the most current solve.
     */
    bool isLOADetected() const { return false; }

    //@}

    //! @name Set methods
    //@{

    //! Set the linear problem that needs to be solved.
    void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; }

    //! Set the parameters the solver manager should use to solve the linear problem.
    void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );

    //@}

    //! @name Reset methods
    //@{
    /*! \brief Performs a reset of the solver manager specified by the \c ResetType.  This informs the
     *  solver manager that the solver should prepare for the next call to solve by resetting certain elements
     *  of the iterative solver strategy.
     */
    void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
    //@}

    //! @name Solver application methods
    //@{

    /*! \brief This method performs possibly repeated calls to the underlying linear solver's
     *         iterate() routine until the problem has been solved (as decided by the solver manager)
     *         or the solver manager decides to quit.
     *
     * This method calls BlockCGIter::iterate() or CGIter::iterate(), which will return either because a
     * specially constructed status test evaluates to ::Passed or an std::exception is thrown.
     *
     * A return from BlockCGIter::iterate() signifies one of the following scenarios:
     *    - the maximum number of iterations has been exceeded. In this scenario, the current solutions
     *      to the linear system will be placed in the linear problem and return ::Unconverged.
     *    - global convergence has been met. In this case, the current solutions to the linear system
     *      will be placed in the linear problem and the solver manager will return ::Converged
     *
     * \returns ::ReturnType specifying:
     *     - ::Converged: the linear problem was solved to the specification required by the solver manager.
     *     - ::Unconverged: the linear problem was not solved to the specification desired by the solver manager.
     */
    ReturnType solve();

    //@}

    /** \name Overridden from Teuchos::Describable */
    //@{

    /** \brief Method to return description of the block CG solver manager */
    std::string description() const;

    //@}

  private:

    //! The linear problem to solve.
    Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;

    //! Output manager, that handles printing of different kinds of messages.
    Teuchos::RCP<OutputManager<ScalarType> > printer_;
    //! Output stream to which the output manager prints.
    Teuchos::RCP<std::ostream> outputStream_;

    /// \brief Aggregate stopping criterion.
    ///
    /// This is an OR combination of the maximum iteration count test
    /// (\c maxIterTest_) and convergence test (\c convTest_).
    Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;

    //! Maximum iteration count stopping criterion.
    Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;

    //! Convergence stopping criterion.
    Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTest_;

    //! Output "status test" that controls all the other status tests.
    Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;

    //! Orthogonalization manager.
    Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;

    //! Current parameter list.
    Teuchos::RCP<Teuchos::ParameterList> params_;

    //
    // Default solver parameters.
    //
    static const MagnitudeType convtol_default_;
    static const MagnitudeType orthoKappa_default_;
    static const int maxIters_default_;
    static const bool adaptiveBlockSize_default_;
    static const bool showMaxResNormOnly_default_;
    static const int blockSize_default_;
    static const int verbosity_default_;
    static const int outputStyle_default_;
    static const int outputFreq_default_;
    static const std::string label_default_;
    static const std::string orthoType_default_;
    static const Teuchos::RCP<std::ostream> outputStream_default_;

    //
    // Current solver parameters and other values.
    //

    //! Convergence tolerance (read from parameter list).
    MagnitudeType convtol_;

    //! Orthogonalization parameter (read from parameter list).
    MagnitudeType orthoKappa_;

    /// \brief Tolerance achieved by the last \c solve() invocation.
    ///
    /// This is the maximum over all right-hand sides' achieved
    /// convergence tolerances, and is set whether or not the solve
    /// actually managed to achieve the desired convergence tolerance.
    MagnitudeType achievedTol_;

    //! Maximum iteration count (read from parameter list).
    int maxIters_;

    //! Number of iterations taken by the last \c solve() invocation.
    int numIters_;

    int blockSize_, verbosity_, outputStyle_, outputFreq_;
    bool adaptiveBlockSize_, showMaxResNormOnly_;
    std::string orthoType_;

    //! Prefix label for all the timers.
    std::string label_;

    //! Solve timer.
    Teuchos::RCP<Teuchos::Time> timerSolve_;

    //! Whether or not the parameters have been set (via \c setParameters()).
    bool isSet_;
  };


// Default solver values.
template<class ScalarType, class MV, class OP>
const typename BlockCGSolMgr<ScalarType,MV,OP,true>::MagnitudeType BlockCGSolMgr<ScalarType,MV,OP,true>::convtol_default_ = 1e-8;

template<class ScalarType, class MV, class OP>
const typename BlockCGSolMgr<ScalarType,MV,OP,true>::MagnitudeType BlockCGSolMgr<ScalarType,MV,OP,true>::orthoKappa_default_ = -1.0;

template<class ScalarType, class MV, class OP>
const int BlockCGSolMgr<ScalarType,MV,OP,true>::maxIters_default_ = 1000;

template<class ScalarType, class MV, class OP>
const bool BlockCGSolMgr<ScalarType,MV,OP,true>::adaptiveBlockSize_default_ = true;

template<class ScalarType, class MV, class OP>
const bool BlockCGSolMgr<ScalarType,MV,OP,true>::showMaxResNormOnly_default_ = false;

template<class ScalarType, class MV, class OP>
const int BlockCGSolMgr<ScalarType,MV,OP,true>::blockSize_default_ = 1;

template<class ScalarType, class MV, class OP>
const int BlockCGSolMgr<ScalarType,MV,OP,true>::verbosity_default_ = Belos::Errors;

template<class ScalarType, class MV, class OP>
const int BlockCGSolMgr<ScalarType,MV,OP,true>::outputStyle_default_ = Belos::General;

template<class ScalarType, class MV, class OP>
const int BlockCGSolMgr<ScalarType,MV,OP,true>::outputFreq_default_ = -1;

template<class ScalarType, class MV, class OP>
const std::string BlockCGSolMgr<ScalarType,MV,OP,true>::label_default_ = "Belos";

template<class ScalarType, class MV, class OP>
const std::string BlockCGSolMgr<ScalarType,MV,OP,true>::orthoType_default_ = "DGKS";

template<class ScalarType, class MV, class OP>
const Teuchos::RCP<std::ostream> BlockCGSolMgr<ScalarType,MV,OP,true>::outputStream_default_ = Teuchos::rcp(&std::cout,false);


// Empty Constructor
template<class ScalarType, class MV, class OP>
BlockCGSolMgr<ScalarType,MV,OP,true>::BlockCGSolMgr() :
  outputStream_(outputStream_default_),
  convtol_(convtol_default_),
  orthoKappa_(orthoKappa_default_),
  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
  maxIters_(maxIters_default_),
  numIters_(0),
  blockSize_(blockSize_default_),
  verbosity_(verbosity_default_),
  outputStyle_(outputStyle_default_),
  outputFreq_(outputFreq_default_),
  adaptiveBlockSize_(adaptiveBlockSize_default_),
  showMaxResNormOnly_(showMaxResNormOnly_default_),
  orthoType_(orthoType_default_),
  label_(label_default_),
  isSet_(false)
{}


// Basic Constructor
template<class ScalarType, class MV, class OP>
BlockCGSolMgr<ScalarType,MV,OP,true>::
BlockCGSolMgr(const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
              const Teuchos::RCP<Teuchos::ParameterList> &pl) :
  problem_(problem),
  outputStream_(outputStream_default_),
  convtol_(convtol_default_),
  orthoKappa_(orthoKappa_default_),
  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
  maxIters_(maxIters_default_),
  numIters_(0),
  blockSize_(blockSize_default_),
  verbosity_(verbosity_default_),
  outputStyle_(outputStyle_default_),
  outputFreq_(outputFreq_default_),
  adaptiveBlockSize_(adaptiveBlockSize_default_),
  showMaxResNormOnly_(showMaxResNormOnly_default_),
  orthoType_(orthoType_default_),
  label_(label_default_),
  isSet_(false)
{
  TEUCHOS_TEST_FOR_EXCEPTION(problem_.is_null(), std::invalid_argument,
    "BlockCGSolMgr's constructor requires a nonnull LinearProblem instance.");

  // If the user passed in a nonnull parameter list, set parameters.
  // Otherwise, the next solve() call will use default parameters,
  // unless the user calls setParameters() first.
  if (! pl.is_null()) {
    setParameters (pl);
  }
}

template<class ScalarType, class MV, class OP>
void
BlockCGSolMgr<ScalarType,MV,OP,true>::
setParameters (const Teuchos::RCP<Teuchos::ParameterList> &params)
{
  // Create the internal parameter list if one doesn't already exist.
  if (params_ == Teuchos::null) {
    params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
  }
  else {
    params->validateParameters(*getValidParameters());
  }

  // Check for maximum number of iterations
  if (params->isParameter("Maximum Iterations")) {
    maxIters_ = params->get("Maximum Iterations",maxIters_default_);

    // Update parameter in our list and in status test.
    params_->set("Maximum Iterations", maxIters_);
    if (maxIterTest_!=Teuchos::null)
      maxIterTest_->setMaxIters( maxIters_ );
  }

  // Check for blocksize
  if (params->isParameter("Block Size")) {
    blockSize_ = params->get("Block Size",blockSize_default_);
    TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
                       "Belos::BlockCGSolMgr: \"Block Size\" must be strictly positive.");

    // Update parameter in our list.
    params_->set("Block Size", blockSize_);
  }

  // Check if the blocksize should be adaptive
  if (params->isParameter("Adaptive Block Size")) {
    adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);

    // Update parameter in our list.
    params_->set("Adaptive Block Size", adaptiveBlockSize_);
  }

  // Check to see if the timer label changed.
  if (params->isParameter("Timer Label")) {
    std::string tempLabel = params->get("Timer Label", label_default_);

    // Update parameter in our list and solver timer
    if (tempLabel != label_) {
      label_ = tempLabel;
      params_->set("Timer Label", label_);
      std::string solveLabel = label_ + ": BlockCGSolMgr total solve time";
#ifdef BELOS_TEUCHOS_TIME_MONITOR
      timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
#endif
      if (ortho_ != Teuchos::null) {
        ortho_->setLabel( label_ );
      }
    }
  }

  // Check if the orthogonalization changed.
  if (params->isParameter("Orthogonalization")) {
    std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
    TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS",
                        std::invalid_argument,
                        "Belos::BlockCGSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
    if (tempOrthoType != orthoType_) {
      orthoType_ = tempOrthoType;
      // Create orthogonalization manager
      if (orthoType_=="DGKS") {
        if (orthoKappa_ <= 0) {
          ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
        }
        else {
          ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
          Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
        }
      }
      else if (orthoType_=="ICGS") {
        ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
      }
      else if (orthoType_=="IMGS") {
        ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
      }
    }
  }

  // Check which orthogonalization constant to use.
  if (params->isParameter("Orthogonalization Constant")) {
    orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);

    // Update parameter in our list.
    params_->set("Orthogonalization Constant",orthoKappa_);
    if (orthoType_=="DGKS") {
      if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
        Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
      }
    }
  }

  // Check for a change in verbosity level
  if (params->isParameter("Verbosity")) {
    if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
      verbosity_ = params->get("Verbosity", verbosity_default_);
    } else {
      verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
    }

    // Update parameter in our list.
    params_->set("Verbosity", verbosity_);
    if (printer_ != Teuchos::null)
      printer_->setVerbosity(verbosity_);
  }

  // Check for a change in output style
  if (params->isParameter("Output Style")) {
    if (Teuchos::isParameterType<int>(*params,"Output Style")) {
      outputStyle_ = params->get("Output Style", outputStyle_default_);
    } else {
      outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
    }

    // Update parameter in our list.
    params_->set("Output Style", outputStyle_);
    outputTest_ = Teuchos::null;
  }

  // output stream
  if (params->isParameter("Output Stream")) {
    outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");

    // Update parameter in our list.
    params_->set("Output Stream", outputStream_);
    if (printer_ != Teuchos::null)
      printer_->setOStream( outputStream_ );
  }

  // frequency level
  if (verbosity_ & Belos::StatusTestDetails) {
    if (params->isParameter("Output Frequency")) {
      outputFreq_ = params->get("Output Frequency", outputFreq_default_);
    }

    // Update parameter in out list and output status test.
    params_->set("Output Frequency", outputFreq_);
    if (outputTest_ != Teuchos::null)
      outputTest_->setOutputFrequency( outputFreq_ );
  }

  // Create output manager if we need to.
  if (printer_ == Teuchos::null) {
    printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
  }

  // Convergence
  typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestResNorm_t;

  // Check for convergence tolerance
  if (params->isParameter("Convergence Tolerance")) {
    convtol_ = params->get("Convergence Tolerance",convtol_default_);

    // Update parameter in our list and residual tests.
    params_->set("Convergence Tolerance", convtol_);
    if (convTest_ != Teuchos::null)
      convTest_->setTolerance( convtol_ );
  }

  if (params->isParameter("Show Maximum Residual Norm Only")) {
    showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");

    // Update parameter in our list and residual tests
    params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
    if (convTest_ != Teuchos::null)
      convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
  }

  // Create status tests if we need to.

  // Basic test checks maximum iterations and native residual.
  if (maxIterTest_ == Teuchos::null)
    maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );

  // Implicit residual test, using the native residual to determine if convergence was achieved.
  if (convTest_ == Teuchos::null)
    convTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_, 1 ) );

  if (sTest_ == Teuchos::null)
    sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );

  if (outputTest_ == Teuchos::null) {

    // Create the status test output class.
    // This class manages and formats the output from the status test.
    StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
    outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );

    // Set the solver string for the output test
    std::string solverDesc = " Block CG ";
    outputTest_->setSolverDesc( solverDesc );

  }

  // Create orthogonalization manager if we need to.
  if (ortho_ == Teuchos::null) {
    if (orthoType_=="DGKS") {
      if (orthoKappa_ <= 0) {
        ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
      }
      else {
        ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
        Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
      }
    }
    else if (orthoType_=="ICGS") {
      ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
    }
    else if (orthoType_=="IMGS") {
      ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
    }
    else {
      TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
                         "Belos::BlockCGSolMgr(): Invalid orthogonalization type.");
    }
  }

  // Create the timer if we need to.
  if (timerSolve_ == Teuchos::null) {
    std::string solveLabel = label_ + ": BlockCGSolMgr total solve time";
#ifdef BELOS_TEUCHOS_TIME_MONITOR
    timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
#endif
  }

  // Inform the solver manager that the current parameters were set.
  isSet_ = true;
}


template<class ScalarType, class MV, class OP>
Teuchos::RCP<const Teuchos::ParameterList>
BlockCGSolMgr<ScalarType,MV,OP,true>::getValidParameters() const
{
  static Teuchos::RCP<const Teuchos::ParameterList> validPL;

  // Set all the valid parameters and their default values.
  if(is_null(validPL)) {
    Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
    pl->set("Convergence Tolerance", convtol_default_,
      "The relative residual tolerance that needs to be achieved by the\n"
      "iterative solver in order for the linear system to be declared converged.");
    pl->set("Maximum Iterations", maxIters_default_,
      "The maximum number of block iterations allowed for each\n"
      "set of RHS solved.");
    pl->set("Block Size", blockSize_default_,
      "The number of vectors in each block.");
    pl->set("Adaptive Block Size", adaptiveBlockSize_default_,
      "Whether the solver manager should adapt to the block size\n"
      "based on the number of RHS to solve.");
    pl->set("Verbosity", verbosity_default_,
      "What type(s) of solver information should be outputted\n"
      "to the output stream.");
    pl->set("Output Style", outputStyle_default_,
      "What style is used for the solver information outputted\n"
      "to the output stream.");
    pl->set("Output Frequency", outputFreq_default_,
      "How often convergence information should be outputted\n"
      "to the output stream.");
    pl->set("Output Stream", outputStream_default_,
      "A reference-counted pointer to the output stream where all\n"
      "solver output is sent.");
    pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
      "When convergence information is printed, only show the maximum\n"
      "relative residual norm when the block size is greater than one.");
    pl->set("Timer Label", label_default_,
      "The string to use as a prefix for the timer labels.");
    //  pl->set("Restart Timers", restartTimers_);
    pl->set("Orthogonalization", orthoType_default_,
      "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
    pl->set("Orthogonalization Constant",orthoKappa_default_,
      "The constant used by DGKS orthogonalization to determine\n"
      "whether another step of classical Gram-Schmidt is necessary.");
    validPL = pl;
  }
  return validPL;
}


// solve()
template<class ScalarType, class MV, class OP>
ReturnType BlockCGSolMgr<ScalarType,MV,OP,true>::solve() {
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Teuchos::rcp_const_cast;
  using Teuchos::rcp_dynamic_cast;

  // Set the current parameters if they were not set before.  NOTE:
  // This may occur if the user generated the solver manager with the
  // default constructor and then didn't set any parameters using
  // setParameters().
  if (!isSet_) {
    setParameters(Teuchos::parameterList(*getValidParameters()));
  }

  Teuchos::BLAS<int,ScalarType> blas;
  Teuchos::LAPACK<int,ScalarType> lapack;

  TEUCHOS_TEST_FOR_EXCEPTION( !problem_->isProblemSet(),
    BlockCGSolMgrLinearProblemFailure,
    "Belos::BlockCGSolMgr::solve(): Linear problem is not ready, setProblem() "
    "has not been called.");

  // Create indices for the linear systems to be solved.
  int startPtr = 0;
  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;

  std::vector<int> currIdx, currIdx2;
  //  If an adaptive block size is allowed then only the linear
  //  systems that need to be solved are solved.  Otherwise, the index
  //  set is generated that informs the linear problem that some
  //  linear systems are augmented.
  if ( adaptiveBlockSize_ ) {
    blockSize_ = numCurrRHS;
    currIdx.resize( numCurrRHS  );
    currIdx2.resize( numCurrRHS  );
    for (int i=0; i<numCurrRHS; ++i)
      { currIdx[i] = startPtr+i; currIdx2[i]=i; }

  }
  else {
    currIdx.resize( blockSize_ );
    currIdx2.resize( blockSize_ );
    for (int i=0; i<numCurrRHS; ++i)
      { currIdx[i] = startPtr+i; currIdx2[i]=i; }
    for (int i=numCurrRHS; i<blockSize_; ++i)
      { currIdx[i] = -1; currIdx2[i] = i; }
  }

  // Inform the linear problem of the current linear system to solve.
  problem_->setLSIndex( currIdx );

  ////////////////////////////////////////////////////////////////////////////
  // Set up the parameter list for the Iteration subclass.
  Teuchos::ParameterList plist;
  plist.set("Block Size",blockSize_);

  // Reset the output status test (controls all the other status tests).
  outputTest_->reset();

  // Assume convergence is achieved, then let any failed convergence
  // set this to false.  "Innocent until proven guilty."
  bool isConverged = true;

  ////////////////////////////////////////////////////////////////////////////
  // Set up the BlockCG Iteration subclass.

  RCP<CGIteration<ScalarType,MV,OP> > block_cg_iter;
  if (blockSize_ == 1) {
    // Standard (nonblock) CG is faster for the special case of a
    // block size of 1.
    block_cg_iter =
      rcp (new CGIter<ScalarType,MV,OP> (problem_, printer_,
                                         outputTest_, plist));
  } else {
    block_cg_iter =
      rcp (new BlockCGIter<ScalarType,MV,OP> (problem_, printer_, outputTest_,
                                              ortho_, plist));
  }


  // Enter solve() iterations
  {
#ifdef BELOS_TEUCHOS_TIME_MONITOR
    Teuchos::TimeMonitor slvtimer(*timerSolve_);
#endif

    while ( numRHS2Solve > 0 ) {
      //
      // Reset the active / converged vectors from this block
      std::vector<int> convRHSIdx;
      std::vector<int> currRHSIdx( currIdx );
      currRHSIdx.resize(numCurrRHS);

      // Reset the number of iterations.
      block_cg_iter->resetNumIters();

      // Reset the number of calls that the status test output knows about.
      outputTest_->resetNumCalls();

      // Get the current residual for this block of linear systems.
      RCP<MV> R_0 = MVT::CloneViewNonConst( *(rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );

      // Set the new state and initialize the solver.
      CGIterationState<ScalarType,MV> newstate;
      newstate.R = R_0;
      block_cg_iter->initializeCG(newstate);

      while(1) {

        // tell block_cg_iter to iterate
        try {
          block_cg_iter->iterate();
          //
          // Check whether any of the linear systems converged.
          //
          if (convTest_->getStatus() == Passed) {
            // At least one of the linear system(s) converged.
            //
            // Get the column indices of the linear systems that converged.
            typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
            std::vector<int> convIdx =
              rcp_dynamic_cast<conv_test_type>(convTest_)->convIndices();

            // If the number of converged linear systems equals the
            // number of linear systems currently being solved, then
            // we are done with this block.
            if (convIdx.size() == currRHSIdx.size())
              break;  // break from while(1){block_cg_iter->iterate()}

            // Inform the linear problem that we are finished with
            // this current linear system.
            problem_->setCurrLS();

            // Reset currRHSIdx to contain the right-hand sides that
            // are left to converge for this block.
            int have = 0;
            std::vector<int> unconvIdx(currRHSIdx.size());
            for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
              bool found = false;
              for (unsigned int j=0; j<convIdx.size(); ++j) {
                if (currRHSIdx[i] == convIdx[j]) {
                  found = true;
                  break;
                }
              }
              if (!found) {
                currIdx2[have] = currIdx2[i];
                currRHSIdx[have++] = currRHSIdx[i];
              }
              else {
              }
            }
            currRHSIdx.resize(have);
            currIdx2.resize(have);

            // Set the remaining indices after deflation.
            problem_->setLSIndex( currRHSIdx );

            // Get the current residual vector.
            std::vector<MagnitudeType> norms;
            R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
            for (int i=0; i<have; ++i) { currIdx2[i] = i; }

            // Set the new blocksize for the solver.
            block_cg_iter->setBlockSize( have );

            // Set the new state and initialize the solver.
            CGIterationState<ScalarType,MV> defstate;
            defstate.R = R_0;
            block_cg_iter->initializeCG(defstate);
          }
          //
          // None of the linear systems converged.  Check whether the
          // maximum iteration count was reached.
          //
          else if (maxIterTest_->getStatus() == Passed) {
            isConverged = false; // None of the linear systems converged.
            break;  // break from while(1){block_cg_iter->iterate()}
          }
          //
          // iterate() returned, but none of our status tests Passed.
          // This indicates a bug.
          //
          else {
            TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
              "Belos::BlockCGSolMgr::solve(): Neither the convergence test nor "
              "the maximum iteration count test passed.  Please report this bug "
              "to the Belos developers.");
          }
        }
        catch (const std::exception &e) {
          std::ostream& err = printer_->stream (Errors);
          err << "Error! Caught std::exception in CGIteration::iterate() at "
              << "iteration " << block_cg_iter->getNumIters() << std::endl
              << e.what() << std::endl;
          throw;
        }
      }

      // Inform the linear problem that we are finished with this
      // block linear system.
      problem_->setCurrLS();

      // Update indices for the linear systems to be solved.
      startPtr += numCurrRHS;
      numRHS2Solve -= numCurrRHS;
      if ( numRHS2Solve > 0 ) {
        numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;

        if ( adaptiveBlockSize_ ) {
          blockSize_ = numCurrRHS;
          currIdx.resize( numCurrRHS  );
          currIdx2.resize( numCurrRHS  );
          for (int i=0; i<numCurrRHS; ++i)
            { currIdx[i] = startPtr+i; currIdx2[i] = i; }
        }
        else {
          currIdx.resize( blockSize_ );
          currIdx2.resize( blockSize_ );
          for (int i=0; i<numCurrRHS; ++i)
            { currIdx[i] = startPtr+i; currIdx2[i] = i; }
          for (int i=numCurrRHS; i<blockSize_; ++i)
            { currIdx[i] = -1; currIdx2[i] = i; }
        }
        // Set the next indices.
        problem_->setLSIndex( currIdx );

        // Set the new blocksize for the solver.
        block_cg_iter->setBlockSize( blockSize_ );
      }
      else {
        currIdx.resize( numRHS2Solve );
      }

    }// while ( numRHS2Solve > 0 )

  }

  // print final summary
  sTest_->print( printer_->stream(FinalSummary) );

  // print timing information
#ifdef BELOS_TEUCHOS_TIME_MONITOR
  // Calling summarize() requires communication in general, so don't
  // call it unless the user wants to print out timing details.
  // summarize() will do all the work even if it's passed a "black
  // hole" output stream.
  if (verbosity_ & TimingDetails) {
    Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
  }
#endif

  // Save the iteration count for this solve.
  numIters_ = maxIterTest_->getNumIters();

  // Save the convergence test value ("achieved tolerance") for this solve.
  {
    typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
    // testValues is nonnull and not persistent.
    const std::vector<MagnitudeType>* pTestValues =
      rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();

    TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
      "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
      "method returned NULL.  Please report this bug to the Belos developers.");

    TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
      "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
      "method returned a vector of length zero.  Please report this bug to the "
      "Belos developers.");

    // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
    // achieved tolerances for all vectors in the current solve(), or
    // just for the vectors from the last deflation?
    achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
  }

  if (!isConverged) {
    return Unconverged; // return from BlockCGSolMgr::solve()
  }
  return Converged; // return from BlockCGSolMgr::solve()
}

//  This method requires the solver manager to return a std::string that describes itself.
template<class ScalarType, class MV, class OP>
std::string BlockCGSolMgr<ScalarType,MV,OP,true>::description() const
{
  std::ostringstream oss;
  oss << "Belos::BlockCGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
  oss << "{";
  oss << "Ortho Type='"<<orthoType_<<"\', Block Size=" << blockSize_;
  oss << "}";
  return oss.str();
}

} // end Belos namespace

#endif /* BELOS_BLOCK_CG_SOLMGR_HPP */