This file is indexed.

/usr/include/trilinos/Thyra_BelosLinearOpWithSolve_def.hpp is in libtrilinos-stratimikos-dev 12.4.2-2.

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
/*
// @HEADER
// ***********************************************************************
//
//         Stratimikos: Thyra-based strategies for linear solvers
//                Copyright (2006) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// 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 Roscoe A. Bartlett (rabartl@sandia.gov)
//
// ***********************************************************************
// @HEADER
*/


#ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
#define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP

#include "Thyra_BelosLinearOpWithSolve_decl.hpp"
#include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp"
#include "Thyra_LinearOpWithSolveHelpers.hpp"
#include "Teuchos_DebugDefaultAsserts.hpp"
#include "Teuchos_Assert.hpp"
#include "Teuchos_TimeMonitor.hpp"
#include "Teuchos_TypeTraits.hpp"

namespace {
  // Set the Belos solver's parameter list to scale its residual norms
  // in the specified way.
  //
  // We break this out in a separate function because the parameters
  // to set depend on which parameters the Belos solver supports.  Not
  // all Belos solvers support both the "Implicit Residual Scaling"
  // and "Explicit Residual Scaling" parameters, so we have to check
  // the solver's list of valid parameters for the existence of these.
  //
  // Scaling options: Belos lets you decide whether the solver will
  // scale residual norms by the (left-)preconditioned initial
  // residual norms (residualScalingType = "Norm of Initial
  // Residual"), or by the unpreconditioned initial residual norms
  // (residualScalingType = "Norm of RHS").  Usually you want to scale
  // by the unpreconditioned initial residual norms.  This is because
  // preconditioning is just an optimization, and you really want to
  // make ||B - AX|| small, rather than ||M B - M (A X)||.  If you're
  // measuring ||B - AX|| and scaling by the initial residual, you
  // should use the unpreconditioned initial residual to match it.
  //
  // Note, however, that the implicit residual test computes
  // left-preconditioned residuals, if a left preconditioner was
  // provided.  That's OK because when Belos solvers (at least the
  // GMRES variants) are given a left preconditioner, they first check
  // the implicit residuals.  If those converge, they then check the
  // explicit residuals.  The explicit residual test does _not_ apply
  // the left preconditioner when computing the residual.  The
  // implicit residual test is just an optimization so that Belos
  // doesn't have to compute explicit residuals B - A*X at every
  // iteration.  This is why we use the same scaling factor for both
  // the implicit and explicit residuals.
  //
  // Arguments:
  //
  // solverParams [in/out] Parameters for the current solve.
  //
  // solverValidParams [in] Valid parameter list for the Belos solver.
  //   Result of calling the solver's getValidParameters() method.
  //
  // residualScalingType [in] String describing how the solver should
  //   scale residuals.  Valid values include "Norm of RHS" and "Norm
  //   of Initial Residual" (these are the only two options this file
  //   currently uses, though Belos offers other options).
  void
  setResidualScalingType (const Teuchos::RCP<Teuchos::ParameterList>& solverParams,
                          const Teuchos::RCP<const Teuchos::ParameterList>& solverValidParams,
                          const std::string& residualScalingType)
  {
    // Many Belos solvers (especially the GMRES variants) define both
    // "Implicit Residual Scaling" and "Explicit Residual Scaling"
    // options.
    //
    // "Implicit" means "the left-preconditioned approximate
    // a.k.a. 'recursive' residual as computed by the Krylov method."
    //
    // "Explicit" means ||B - A*X||, the unpreconditioned, "exact"
    // residual.
    //
    // Belos' GMRES implementations chain these two tests in sequence.
    // Implicit comes first, and explicit is not evaluated unless
    // implicit passes.  In some cases (e.g., no left preconditioner),
    // GMRES _only_ uses the implicit tests.  This means that only
    // setting "Explicit Residual Scaling" won't change the solver's
    // behavior.  Stratimikos tends to prefer using a right
    // preconditioner, in which case setting only the "Explicit
    // Residual Scaling" argument has no effect.  Furthermore, if
    // "Explicit Residual Scaling" is set to something other than the
    // default (initial residual norm), without "Implicit Residual
    // Scaling" getting the same setting, then the implicit residual
    // test will be using a radically different scaling factor than
    // the user wanted.
    //
    // Not all Belos solvers support both options.  We check the
    // solver's valid parameter list first before attempting to set
    // the option.
    if (solverValidParams->isParameter ("Implicit Residual Scaling")) {
      solverParams->set ("Implicit Residual Scaling", residualScalingType);
    }
    if (solverValidParams->isParameter ("Explicit Residual Scaling")) {
      solverParams->set ("Explicit Residual Scaling", residualScalingType);
    }
  }

} // namespace (anonymous)


namespace Thyra {


// Constructors/initializers/accessors


template<class Scalar>
BelosLinearOpWithSolve<Scalar>::BelosLinearOpWithSolve()
  :convergenceTestFrequency_(-1),
  isExternalPrec_(false),
  supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED),
  defaultTol_ (-1.0)
{}


template<class Scalar>
void BelosLinearOpWithSolve<Scalar>::initialize(
  const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp,
  const RCP<Teuchos::ParameterList> &solverPL,
  const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
  const RCP<const PreconditionerBase<Scalar> > &prec,
  const bool isExternalPrec_in,
  const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
  const ESupportSolveUse &supportSolveUse_in,
  const int convergenceTestFrequency
  )
{
  using Teuchos::as;
  using Teuchos::TypeNameTraits;
  using Teuchos::Exceptions::InvalidParameterType;
  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;

  this->setLinePrefix("BELOS/T");
  // ToDo: Validate input
  lp_ = lp;
  solverPL_ = solverPL;
  iterativeSolver_ = iterativeSolver;
  fwdOpSrc_ = fwdOpSrc;
  prec_ = prec;
  isExternalPrec_ = isExternalPrec_in;
  approxFwdOpSrc_ = approxFwdOpSrc;
  supportSolveUse_ = supportSolveUse_in;
  convergenceTestFrequency_ = convergenceTestFrequency;
  // Check if "Convergence Tolerance" is in the solver parameter list.  If
  // not, use the default from the solver.
  if ( !is_null(solverPL_) ) {
    if (solverPL_->isParameter("Convergence Tolerance")) {

      // Stratimikos prefers tolerances as double, no matter the
      // Scalar type.  However, we also want it to accept the
      // tolerance as magnitude_type, for example float if the Scalar
      // type is float or std::complex<float>.
      if (solverPL_->isType<double> ("Convergence Tolerance")) {
        defaultTol_ =
          as<magnitude_type> (solverPL_->get<double> ("Convergence Tolerance"));
      }
      else if (Teuchos::TypeTraits::is_same<double, magnitude_type>::value) {
        // magnitude_type == double in this case, and we've already
        // checked double above.
        TEUCHOS_TEST_FOR_EXCEPTION(
          true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
          "The \"Convergence Tolerance\" parameter, which you provided, must "
          "have type double (the type of the magnitude of Scalar = double).");
      }
      else if (solverPL_->isType<magnitude_type> ("Convergence Tolerance")) {
        defaultTol_ = solverPL_->get<magnitude_type> ("Convergence Tolerance");
      }
      else {
        // Throwing InvalidParameterType ensures that the exception's
        // type is consistent both with what this method would have
        // thrown before for an unrecognized type, and with what the
        // user expects in general when the parameter doesn't have the
        // right type.
        TEUCHOS_TEST_FOR_EXCEPTION(
          true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
          "The \"Convergence Tolerance\" parameter, which you provided, must "
          "have type double (preferred) or the type of the magnitude of Scalar "
          "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
          TypeNameTraits<magnitude_type>::name () << " in this case.  You can "
          "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
      }
    }
  }
  else {
    RCP<const Teuchos::ParameterList> defaultPL =
      iterativeSolver->getValidParameters();

    // Stratimikos prefers tolerances as double, no matter the
    // Scalar type.  However, we also want it to accept the
    // tolerance as magnitude_type, for example float if the Scalar
    // type is float or std::complex<float>.
    if (defaultPL->isType<double> ("Convergence Tolerance")) {
      defaultTol_ =
        as<magnitude_type> (defaultPL->get<double> ("Convergence Tolerance"));
    }
    else if (Teuchos::TypeTraits::is_same<double, magnitude_type>::value) {
      // magnitude_type == double in this case, and we've already
      // checked double above.
      TEUCHOS_TEST_FOR_EXCEPTION(
        true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
        "The \"Convergence Tolerance\" parameter, which you provided, must "
        "have type double (the type of the magnitude of Scalar = double).");
    }
    else if (defaultPL->isType<magnitude_type> ("Convergence Tolerance")) {
      defaultTol_ = defaultPL->get<magnitude_type> ("Convergence Tolerance");
    }
    else {
      // Throwing InvalidParameterType ensures that the exception's
      // type is consistent both with what this method would have
      // thrown before for an unrecognized type, and with what the
      // user expects in general when the parameter doesn't have the
      // right type.
      TEUCHOS_TEST_FOR_EXCEPTION(
        true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
        "The \"Convergence Tolerance\" parameter, which you provided, must "
        "have type double (preferred) or the type of the magnitude of Scalar "
        "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
        TypeNameTraits<magnitude_type>::name () << " in this case.  You can "
        "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
    }
  }
}


template<class Scalar>
RCP<const LinearOpSourceBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::extract_fwdOpSrc()
{
  RCP<const LinearOpSourceBase<Scalar> >
    _fwdOpSrc = fwdOpSrc_;
  fwdOpSrc_ = Teuchos::null;
  return _fwdOpSrc;
}


template<class Scalar>
RCP<const PreconditionerBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::extract_prec()
{
  RCP<const PreconditionerBase<Scalar> >
    _prec = prec_;
  prec_ = Teuchos::null;
  return _prec;
}


template<class Scalar>
bool BelosLinearOpWithSolve<Scalar>::isExternalPrec() const
{
  return isExternalPrec_;
}


template<class Scalar>
RCP<const LinearOpSourceBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::extract_approxFwdOpSrc()
{
  RCP<const LinearOpSourceBase<Scalar> >
    _approxFwdOpSrc = approxFwdOpSrc_;
  approxFwdOpSrc_ = Teuchos::null;
  return _approxFwdOpSrc;
}


template<class Scalar>
ESupportSolveUse BelosLinearOpWithSolve<Scalar>::supportSolveUse() const
{
  return supportSolveUse_;
}


template<class Scalar>
void BelosLinearOpWithSolve<Scalar>::uninitialize(
  RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp,
  RCP<Teuchos::ParameterList> *solverPL,
  RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
  RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
  RCP<const PreconditionerBase<Scalar> > *prec,
  bool *isExternalPrec_in,
  RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
  ESupportSolveUse *supportSolveUse_in
  )
{
  if (lp) *lp = lp_;
  if (solverPL) *solverPL = solverPL_;
  if (iterativeSolver) *iterativeSolver = iterativeSolver_;
  if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
  if (prec) *prec = prec_;
  if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
  if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
  if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;

  lp_ = Teuchos::null;
  solverPL_ = Teuchos::null;
  iterativeSolver_ = Teuchos::null;
  fwdOpSrc_ = Teuchos::null;
  prec_ = Teuchos::null;
  isExternalPrec_ = false;
  approxFwdOpSrc_ = Teuchos::null;
  supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
}


// Overridden from LinearOpBase


template<class Scalar>
RCP< const VectorSpaceBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::range() const
{
  if (!is_null(lp_))
    return lp_->getOperator()->range();
  return Teuchos::null;
}


template<class Scalar>
RCP< const VectorSpaceBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::domain() const
{
  if (!is_null(lp_))
    return lp_->getOperator()->domain();
  return Teuchos::null;
}


template<class Scalar>
RCP<const LinearOpBase<Scalar> >
BelosLinearOpWithSolve<Scalar>::clone() const
{
  return Teuchos::null; // Not supported yet but could be
}


// Overridden from Teuchos::Describable


template<class Scalar>
std::string BelosLinearOpWithSolve<Scalar>::description() const
{
  std::ostringstream oss;
  oss << Teuchos::Describable::description();
  if ( !is_null(lp_) && !is_null(lp_->getOperator()) ) {
    oss << "{";
    oss << "iterativeSolver=\'"<<iterativeSolver_->description()<<"\'";
    oss << ",fwdOp=\'"<<lp_->getOperator()->description()<<"\'";
    if (lp_->getLeftPrec().get())
      oss << ",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<"\'";
    if (lp_->getRightPrec().get())
      oss << ",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<"\'";
    oss << "}";
  }
  // ToDo: Make Belos::SolverManager derive from Teuchos::Describable so
  // that we can get better information.
  return oss.str();
}


template<class Scalar>
void BelosLinearOpWithSolve<Scalar>::describe(
  Teuchos::FancyOStream &out_arg,
  const Teuchos::EVerbosityLevel verbLevel
  ) const
{
  using Teuchos::FancyOStream;
  using Teuchos::OSTab;
  using Teuchos::describe;
  RCP<FancyOStream> out = rcp(&out_arg,false);
  OSTab tab(out);
  switch (verbLevel) {
    case Teuchos::VERB_LOW:
      break;
    case Teuchos::VERB_DEFAULT:
    case Teuchos::VERB_MEDIUM:
      *out << this->description() << std::endl;
      break;
    case Teuchos::VERB_HIGH:
    case Teuchos::VERB_EXTREME:
    {
      *out
        << Teuchos::Describable::description()<< "{"
        << "rangeDim=" << this->range()->dim()
        << ",domainDim=" << this->domain()->dim() << "}\n";
      if (lp_->getOperator().get()) {
        OSTab tab1(out);
        *out
          << "iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
          << "fwdOp = " << describe(*lp_->getOperator(),verbLevel);
        if (lp_->getLeftPrec().get())
          *out << "leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
        if (lp_->getRightPrec().get())
          *out << "rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
      }
      break;
    }
    default:
      TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
  }
}


// protected


// Overridden from LinearOpBase


template<class Scalar>
bool BelosLinearOpWithSolve<Scalar>::opSupportedImpl(EOpTransp M_trans) const
{
  return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
}


template<class Scalar>
void BelosLinearOpWithSolve<Scalar>::applyImpl(
  const EOpTransp M_trans,
  const MultiVectorBase<Scalar> &X,
  const Ptr<MultiVectorBase<Scalar> > &Y,
  const Scalar alpha,
  const Scalar beta
  ) const
{
  ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
}


// Overridden from LinearOpWithSolveBase


template<class Scalar>
bool
BelosLinearOpWithSolve<Scalar>::solveSupportsImpl(EOpTransp M_trans) const
{
  return solveSupportsNewImpl(M_trans, Teuchos::null);
}


template<class Scalar>
bool
BelosLinearOpWithSolve<Scalar>::solveSupportsNewImpl(EOpTransp transp,
  const Ptr<const SolveCriteria<Scalar> > solveCriteria) const
{
  // Only support forward solve right now!
  if (real_trans(transp)==NOTRANS) return true;
  return false; // ToDo: Support adjoint solves!
  // Otherwise, Thyra/Belos now supports every solve criteria type that exists
  // because of the class Thyra::GeneralSolveCriteriaBelosStatusTest!
  return true;
/*
  if (real_trans(M_trans)==NOTRANS) {
    return (
      solveMeasureType.useDefault()
      ||
      solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_RHS)
      ||
      solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_INIT_RESIDUAL)
      );
  }
*/
}


template<class Scalar>
bool
BelosLinearOpWithSolve<Scalar>::solveSupportsSolveMeasureTypeImpl(
  EOpTransp M_trans, const SolveMeasureType& solveMeasureType) const
{
  SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
  return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
}


template<class Scalar>
SolveStatus<Scalar>
BelosLinearOpWithSolve<Scalar>::solveImpl(
  const EOpTransp M_trans,
  const MultiVectorBase<Scalar> &B,
  const Ptr<MultiVectorBase<Scalar> > &X,
  const Ptr<const SolveCriteria<Scalar> > solveCriteria
  ) const
{

  THYRA_FUNC_TIME_MONITOR("Stratimikos: BelosLOWS");

  using Teuchos::rcp;
  using Teuchos::rcpFromRef;
  using Teuchos::rcpFromPtr;
  using Teuchos::FancyOStream;
  using Teuchos::OSTab;
  using Teuchos::ParameterList;
  using Teuchos::parameterList;
  using Teuchos::describe;
  typedef Teuchos::ScalarTraits<Scalar> ST;
  typedef typename ST::magnitudeType ScalarMag;
  Teuchos::Time totalTimer(""), timer("");
  totalTimer.start(true);

  assertSolveSupports(*this, M_trans, solveCriteria);
  // 2010/08/22: rabartl: Bug 4915 ToDo: Move the above into the NIV function
  // solve(...).

  const RCP<FancyOStream> out = this->getOStream();
  const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
  OSTab tab = this->getOSTab();
  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW)) {
    *out << "\nStarting iterations with Belos:\n";
    OSTab tab2(out);
    *out << "Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
    *out << "Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
    *out << "With #Eqns="<<B.range()->dim()<<", #RHSs="<<B.domain()->dim()<<" ...\n";
  }

  //
  // Set RHS and LHS
  //

  bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
  TEUCHOS_TEST_FOR_EXCEPTION(
    ret == false, CatastrophicSolveFailure
    ,"Error, the Belos::LinearProblem could not be set for the current solve!"
    );

  //
  // Set the solution criteria
  //

  // Parameter list for the current solve.
  const RCP<ParameterList> tmpPL = Teuchos::parameterList();

  // The solver's valid parameter list.
  RCP<const ParameterList> validPL = iterativeSolver_->getValidParameters();

  SolveMeasureType solveMeasureType;
  RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> > generalSolveCriteriaBelosStatusTest;
  if (nonnull(solveCriteria)) {
    solveMeasureType = solveCriteria->solveMeasureType;
    const ScalarMag requestedTol = solveCriteria->requestedTol;
    if (solveMeasureType.useDefault()) {
      tmpPL->set("Convergence Tolerance", defaultTol_);
    }
    else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
      if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
        tmpPL->set("Convergence Tolerance", requestedTol);
      }
      else {
        tmpPL->set("Convergence Tolerance", defaultTol_);
      }
      setResidualScalingType (tmpPL, validPL, "Norm of RHS");
    }
    else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
      if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
        tmpPL->set("Convergence Tolerance", requestedTol);
      }
      else {
        tmpPL->set("Convergence Tolerance", defaultTol_);
      }
      setResidualScalingType (tmpPL, validPL, "Norm of Initial Residual");
    }
    else {
      // Set the most generic (and inefficient) solve criteria
      generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
        *solveCriteria, convergenceTestFrequency_);
      // Set the verbosity level (one level down)
      generalSolveCriteriaBelosStatusTest->setOStream(out);
      generalSolveCriteriaBelosStatusTest->setVerbLevel(incrVerbLevel(verbLevel, -1));
      // Set the default convergence tolerance to always converged to allow
      // the above status test to control things.
      tmpPL->set("Convergence Tolerance", 1.0);
    }
    // maximum iterations
    if (nonnull(solveCriteria->extraParameters)) {
      if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,"Maximum Iterations")) {
        tmpPL->set("Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,"Maximum Iterations"));
      }
    }
    // If a preconditioner is on the left, then the implicit residual test
    // scaling should be the preconditioned initial residual.
    if (Teuchos::nonnull(lp_->getLeftPrec()) &&
        validPL->isParameter ("Implicit Residual Scaling"))
      tmpPL->set("Implicit Residual Scaling",
                 "Norm of Preconditioned Initial Residual");
  }
  else {
    // No solveCriteria was even passed in!
    tmpPL->set("Convergence Tolerance", defaultTol_);
  }

  //
  // Solve the linear system
  //

  Belos::ReturnType belosSolveStatus;
  {
    // Write detailed convergence information if requested for levels >= VERB_LOW
    RCP<std::ostream>
      outUsed =
      ( static_cast<int>(verbLevel) >= static_cast<int>(Teuchos::VERB_LOW)
        ? out
        : rcp(new FancyOStream(rcp(new Teuchos::oblackholestream())))
        );
    Teuchos::OSTab tab1(outUsed,1,"BELOS");
    tmpPL->set("Output Stream", outUsed);
    iterativeSolver_->setParameters(tmpPL);
    if (nonnull(generalSolveCriteriaBelosStatusTest)) {
      iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
    }
    try {
      belosSolveStatus = iterativeSolver_->solve();
    }
    catch (Belos::BelosError&) {
      belosSolveStatus = Belos::Unconverged;
    }
  }

  //
  // Report the solve status
  //

  totalTimer.stop();

  SolveStatus<Scalar> solveStatus;

  switch (belosSolveStatus) {
    case Belos::Unconverged: {
      solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
      // Set achievedTol even if the solver did not converge.  This is
      // helpful for things like nonlinear solvers, which might be
      // able to use a partially converged result, and which would
      // like to know the achieved convergence tolerance for use in
      // computing bounds.  It's also helpful for estimating whether a
      // small increase in the maximum iteration count might be
      // helpful next time.
      try {
        // Some solvers might not have implemented achievedTol().
        // The default implementation throws std::runtime_error.
        solveStatus.achievedTol = iterativeSolver_->achievedTol();
      } catch (std::runtime_error&) {
        // Do nothing; use the default value of achievedTol.
      }
      break;
    }
    case Belos::Converged: {
      solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
      if (nonnull(generalSolveCriteriaBelosStatusTest)) {
        // The user set a custom status test.  This means that we
        // should ask the custom status test itself, rather than the
        // Belos solver, what the final achieved convergence tolerance
        // was.
        const ArrayView<const ScalarMag> achievedTol =
          generalSolveCriteriaBelosStatusTest->achievedTol();
        solveStatus.achievedTol = ST::zero();
        for (Ordinal i = 0; i < achievedTol.size(); ++i) {
          solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
        }
      }
      else {
        try {
          // Some solvers might not have implemented achievedTol().
          // The default implementation throws std::runtime_error.
          solveStatus.achievedTol = iterativeSolver_->achievedTol();
        } catch (std::runtime_error&) {
          // Use the default convergence tolerance.  This is a correct
          // upper bound, since we did actually converge.
          solveStatus.achievedTol = tmpPL->get("Convergence Tolerance", defaultTol_);
        }
      }
      break;
    }
    TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT();
  }

  std::ostringstream ossmessage;
  ossmessage
    << "The Belos solver of type \""<<iterativeSolver_->description()
    <<"\" returned a solve status of \""<< toString(solveStatus.solveStatus) << "\""
    << " in " << iterativeSolver_->getNumIters() << " iterations"
    << " with total CPU time of " << totalTimer.totalElapsedTime() << " sec" ;
  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
    *out << "\n" << ossmessage.str() << "\n";

  solveStatus.message = ossmessage.str();

  // Dump the getNumIters() and the achieved convergence tolerance
  // into solveStatus.extraParameters, as the "Belos/Iteration Count"
  // resp. "Belos/Achieved Tolerance" parameters.
  if (solveStatus.extraParameters.is_null()) {
    solveStatus.extraParameters = parameterList ();
  }
  solveStatus.extraParameters->set ("Belos/Iteration Count",
                                    iterativeSolver_->getNumIters());\
  // package independent version of the same
  solveStatus.extraParameters->set ("Iteration Count",
                                    iterativeSolver_->getNumIters());\
  // NOTE (mfh 13 Dec 2011) Though the most commonly used Belos
  // solvers do implement achievedTol(), some Belos solvers currently
  // do not.  In the latter case, if the solver did not converge, the
  // reported achievedTol() value may just be the default "invalid"
  // value -1, and if the solver did converge, the reported value will
  // just be the convergence tolerance (a correct upper bound).
  solveStatus.extraParameters->set ("Belos/Achieved Tolerance",
                                    solveStatus.achievedTol);

//  This information is in the previous line, which is printed anytime the verbosity
//  is not set to Teuchos::VERB_NONE, so I'm commenting this out for now.
//  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE))
//    *out << "\nTotal solve time in Belos = "<<totalTimer.totalElapsedTime()<<" sec\n";

  return solveStatus;

}


}       // end namespace Thyra


#endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP